File indexing completed on 2024-10-16 05:06:28
0001 #include <string>
0002 #include <memory>
0003
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/Utilities/interface/StreamID.h"
0011
0012 #include "DataFormats/Common/interface/View.h"
0013 #include "DataFormats/Common/interface/ValueMap.h"
0014 #include "DataFormats/PatCandidates/interface/Jet.h"
0015 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0016 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0017 #include "DataFormats/NanoAOD/interface/FlatTable.h"
0018
0019 template <typename T>
0020 class JetPFConstituentVarProducer : public edm::stream::EDProducer<> {
0021 public:
0022 explicit JetPFConstituentVarProducer(const edm::ParameterSet&);
0023 ~JetPFConstituentVarProducer() override;
0024 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0025
0026 private:
0027 void produce(edm::Event&, const edm::EventSetup&) override;
0028 void PutValueMapInEvent(edm::Event&, const edm::Handle<edm::View<T>>&, const std::vector<float>&, std::string);
0029
0030 edm::EDGetTokenT<edm::View<T>> jet_token_;
0031 edm::EDGetTokenT<edm::ValueMap<float>> puppi_value_map_token_;
0032
0033 const bool fallback_puppi_weight_;
0034 bool use_puppi_value_map_;
0035 };
0036
0037
0038
0039
0040 template <typename T>
0041 JetPFConstituentVarProducer<T>::JetPFConstituentVarProducer(const edm::ParameterSet& iConfig)
0042 : jet_token_(consumes<edm::View<T>>(iConfig.getParameter<edm::InputTag>("jets"))),
0043 fallback_puppi_weight_(iConfig.getParameter<bool>("fallback_puppi_weight")),
0044 use_puppi_value_map_(false) {
0045
0046
0047
0048 const auto& puppi_value_map_tag = iConfig.getParameter<edm::InputTag>("puppi_value_map");
0049 if (!puppi_value_map_tag.label().empty()) {
0050 puppi_value_map_token_ = consumes<edm::ValueMap<float>>(puppi_value_map_tag);
0051 use_puppi_value_map_ = true;
0052 }
0053 produces<edm::ValueMap<float>>("leadConstNeHadEF");
0054 produces<edm::ValueMap<float>>("leadConstChHadEF");
0055 produces<edm::ValueMap<float>>("leadConstPhotonEF");
0056 produces<edm::ValueMap<float>>("leadConstElectronEF");
0057 produces<edm::ValueMap<float>>("leadConstMuonEF");
0058 produces<edm::ValueMap<float>>("leadConstHFHADEF");
0059 produces<edm::ValueMap<float>>("leadConstHFEMEF");
0060 produces<edm::ValueMap<float>>("leadConstNeHadPuppiWeight");
0061 produces<edm::ValueMap<float>>("leadConstChHadPuppiWeight");
0062 produces<edm::ValueMap<float>>("leadConstPhotonPuppiWeight");
0063 produces<edm::ValueMap<float>>("leadConstElectronPuppiWeight");
0064 produces<edm::ValueMap<float>>("leadConstMuonPuppiWeight");
0065 produces<edm::ValueMap<float>>("leadConstHFHADPuppiWeight");
0066 produces<edm::ValueMap<float>>("leadConstHFEMPuppiWeight");
0067 }
0068
0069 template <typename T>
0070 JetPFConstituentVarProducer<T>::~JetPFConstituentVarProducer() {}
0071
0072 template <typename T>
0073 void JetPFConstituentVarProducer<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0074
0075 auto jets = iEvent.getHandle(jet_token_);
0076
0077
0078 edm::Handle<edm::ValueMap<float>> puppi_value_map_;
0079 if (use_puppi_value_map_) {
0080 iEvent.getByToken(puppi_value_map_token_, puppi_value_map_);
0081 }
0082
0083 std::vector<float> jet_leadConstNeHadEF(jets->size(), 0.f);
0084 std::vector<float> jet_leadConstNeHadPuppiWeight(jets->size(), 0.f);
0085
0086 std::vector<float> jet_leadConstChHadEF(jets->size(), 0.f);
0087 std::vector<float> jet_leadConstChHadPuppiWeight(jets->size(), 0.f);
0088
0089 std::vector<float> jet_leadConstPhotonEF(jets->size(), 0.f);
0090 std::vector<float> jet_leadConstPhotonPuppiWeight(jets->size(), 0.f);
0091
0092 std::vector<float> jet_leadConstElectronEF(jets->size(), 0.f);
0093 std::vector<float> jet_leadConstElectronPuppiWeight(jets->size(), 0.f);
0094
0095 std::vector<float> jet_leadConstMuonEF(jets->size(), 0.f);
0096 std::vector<float> jet_leadConstMuonPuppiWeight(jets->size(), 0.f);
0097
0098 std::vector<float> jet_leadConstHFHADEF(jets->size(), 0.f);
0099 std::vector<float> jet_leadConstHFHADPuppiWeight(jets->size(), 0.f);
0100
0101 std::vector<float> jet_leadConstHFEMEF(jets->size(), 0.f);
0102 std::vector<float> jet_leadConstHFEMPuppiWeight(jets->size(), 0.f);
0103
0104
0105 for (std::size_t jet_idx = 0; jet_idx < jets->size(); jet_idx++) {
0106 const auto& jet = (*jets)[jet_idx];
0107
0108 float jet_energy_raw = jet.energy();
0109 if constexpr (std::is_same<T, pat::Jet>::value) {
0110 jet_energy_raw = jet.correctedJet(0).energy();
0111 }
0112
0113 reco::CandidatePtr leadConstNeHad;
0114 reco::CandidatePtr leadConstChHad;
0115 reco::CandidatePtr leadConstPhoton;
0116 reco::CandidatePtr leadConstElectron;
0117 reco::CandidatePtr leadConstMuon;
0118 reco::CandidatePtr leadConstHFHAD;
0119 reco::CandidatePtr leadConstHFEM;
0120
0121 float leadConstNeHadPuppiWeight = 1.f;
0122 float leadConstChHadPuppiWeight = 1.f;
0123 float leadConstPhotonPuppiWeight = 1.f;
0124 float leadConstElectronPuppiWeight = 1.f;
0125 float leadConstMuonPuppiWeight = 1.f;
0126 float leadConstHFHADPuppiWeight = 1.f;
0127 float leadConstHFEMPuppiWeight = 1.f;
0128
0129
0130
0131
0132 for (const reco::CandidatePtr& dau : jet.daughterPtrVector()) {
0133 float puppiw = 1.f;
0134
0135
0136
0137
0138 if (use_puppi_value_map_) {
0139 puppiw = (*puppi_value_map_)[dau];
0140 } else if (!fallback_puppi_weight_) {
0141 throw edm::Exception(edm::errors::InvalidReference, "PUPPI value map missing")
0142 << "use fallback_puppi_weight option to use " << puppiw << " for cand as default";
0143 }
0144
0145
0146
0147
0148 if (abs(dau->pdgId()) == 130) {
0149 if (leadConstNeHad.isNull() ||
0150 (puppiw * dau->energy() > leadConstNeHadPuppiWeight * leadConstNeHad->energy())) {
0151 leadConstNeHad = dau;
0152 leadConstNeHadPuppiWeight = puppiw;
0153 }
0154 } else if (abs(dau->pdgId()) == 211) {
0155 if (leadConstChHad.isNull() ||
0156 (puppiw * dau->energy() > leadConstChHadPuppiWeight * leadConstChHad->energy())) {
0157 leadConstChHad = dau;
0158 leadConstChHadPuppiWeight = puppiw;
0159 }
0160 } else if (abs(dau->pdgId()) == 22) {
0161 if (leadConstPhoton.isNull() ||
0162 (puppiw * dau->energy() > leadConstPhotonPuppiWeight * leadConstPhoton->energy())) {
0163 leadConstPhoton = dau;
0164 leadConstPhotonPuppiWeight = puppiw;
0165 }
0166 } else if (abs(dau->pdgId()) == 11) {
0167 if (leadConstElectron.isNull() ||
0168 (puppiw * dau->energy() > leadConstElectronPuppiWeight * leadConstElectron->energy())) {
0169 leadConstElectron = dau;
0170 leadConstElectronPuppiWeight = puppiw;
0171 }
0172 } else if (abs(dau->pdgId()) == 13) {
0173 if (leadConstMuon.isNull() || (puppiw * dau->energy() > leadConstMuonPuppiWeight * leadConstMuon->energy())) {
0174 leadConstMuon = dau;
0175 leadConstMuonPuppiWeight = puppiw;
0176 }
0177 } else if (abs(dau->pdgId()) == 1) {
0178 if (leadConstHFHAD.isNull() ||
0179 (puppiw * dau->energy() > leadConstHFHADPuppiWeight * leadConstHFHAD->energy())) {
0180 leadConstHFHAD = dau;
0181 leadConstHFHADPuppiWeight = puppiw;
0182 }
0183 } else if (abs(dau->pdgId()) == 2) {
0184 if (leadConstHFEM.isNull() || (puppiw * dau->energy() > leadConstHFEMPuppiWeight * leadConstHFEM->energy())) {
0185 leadConstHFEM = dau;
0186 leadConstHFEMPuppiWeight = puppiw;
0187 }
0188 }
0189 }
0190
0191 if (leadConstNeHad.isNonnull()) {
0192 jet_leadConstNeHadEF[jet_idx] = (leadConstNeHad->energy() * leadConstNeHadPuppiWeight) / jet_energy_raw;
0193 jet_leadConstNeHadPuppiWeight[jet_idx] = leadConstNeHadPuppiWeight;
0194 }
0195 if (leadConstChHad.isNonnull()) {
0196 jet_leadConstChHadEF[jet_idx] = (leadConstChHad->energy() * leadConstChHadPuppiWeight) / jet_energy_raw;
0197 jet_leadConstChHadPuppiWeight[jet_idx] = leadConstChHadPuppiWeight;
0198 }
0199 if (leadConstPhoton.isNonnull()) {
0200 jet_leadConstPhotonEF[jet_idx] = (leadConstPhoton->energy() * leadConstPhotonPuppiWeight) / jet_energy_raw;
0201 jet_leadConstPhotonPuppiWeight[jet_idx] = leadConstPhotonPuppiWeight;
0202 }
0203 if (leadConstElectron.isNonnull()) {
0204 jet_leadConstElectronEF[jet_idx] = (leadConstElectron->energy() * leadConstElectronPuppiWeight) / jet_energy_raw;
0205 jet_leadConstElectronPuppiWeight[jet_idx] = leadConstElectronPuppiWeight;
0206 }
0207 if (leadConstMuon.isNonnull()) {
0208 jet_leadConstMuonEF[jet_idx] = (leadConstMuon->energy() * leadConstMuonPuppiWeight) / jet_energy_raw;
0209 jet_leadConstMuonPuppiWeight[jet_idx] = leadConstMuonPuppiWeight;
0210 }
0211 if (leadConstHFHAD.isNonnull()) {
0212 jet_leadConstHFHADEF[jet_idx] = (leadConstHFHAD->energy() * leadConstHFHADPuppiWeight) / jet_energy_raw;
0213 jet_leadConstHFHADPuppiWeight[jet_idx] = leadConstHFHADPuppiWeight;
0214 }
0215 if (leadConstHFEM.isNonnull()) {
0216 jet_leadConstHFEMEF[jet_idx] = (leadConstHFEM->energy() * leadConstHFEMPuppiWeight) / jet_energy_raw;
0217 jet_leadConstHFEMPuppiWeight[jet_idx] = leadConstHFEMPuppiWeight;
0218 }
0219 }
0220
0221 PutValueMapInEvent(iEvent, jets, jet_leadConstNeHadEF, "leadConstNeHadEF");
0222 PutValueMapInEvent(iEvent, jets, jet_leadConstNeHadPuppiWeight, "leadConstNeHadPuppiWeight");
0223
0224 PutValueMapInEvent(iEvent, jets, jet_leadConstChHadEF, "leadConstChHadEF");
0225 PutValueMapInEvent(iEvent, jets, jet_leadConstChHadPuppiWeight, "leadConstChHadPuppiWeight");
0226
0227 PutValueMapInEvent(iEvent, jets, jet_leadConstPhotonEF, "leadConstPhotonEF");
0228 PutValueMapInEvent(iEvent, jets, jet_leadConstPhotonPuppiWeight, "leadConstPhotonPuppiWeight");
0229
0230 PutValueMapInEvent(iEvent, jets, jet_leadConstElectronEF, "leadConstElectronEF");
0231 PutValueMapInEvent(iEvent, jets, jet_leadConstElectronPuppiWeight, "leadConstElectronPuppiWeight");
0232
0233 PutValueMapInEvent(iEvent, jets, jet_leadConstMuonEF, "leadConstMuonEF");
0234 PutValueMapInEvent(iEvent, jets, jet_leadConstMuonPuppiWeight, "leadConstMuonPuppiWeight");
0235
0236 PutValueMapInEvent(iEvent, jets, jet_leadConstHFHADEF, "leadConstHFHADEF");
0237 PutValueMapInEvent(iEvent, jets, jet_leadConstHFHADPuppiWeight, "leadConstHFHADPuppiWeight");
0238
0239 PutValueMapInEvent(iEvent, jets, jet_leadConstHFEMEF, "leadConstHFEMEF");
0240 PutValueMapInEvent(iEvent, jets, jet_leadConstHFEMPuppiWeight, "leadConstHFEMPuppiWeight");
0241 }
0242 template <typename T>
0243 void JetPFConstituentVarProducer<T>::PutValueMapInEvent(edm::Event& iEvent,
0244 const edm::Handle<edm::View<T>>& coll,
0245 const std::vector<float>& vec_var,
0246 std::string VMName) {
0247 std::unique_ptr<edm::ValueMap<float>> VM(new edm::ValueMap<float>());
0248 edm::ValueMap<float>::Filler fillerVM(*VM);
0249 fillerVM.insert(coll, vec_var.begin(), vec_var.end());
0250 fillerVM.fill();
0251 iEvent.put(std::move(VM), VMName);
0252 }
0253
0254 template <typename T>
0255 void JetPFConstituentVarProducer<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0256 edm::ParameterSetDescription desc;
0257 desc.add<edm::InputTag>("jets", edm::InputTag("finalJetsPuppi"));
0258 desc.add<edm::InputTag>("puppi_value_map", edm::InputTag("packedpuppi"));
0259 desc.add<bool>("fallback_puppi_weight", false);
0260 descriptions.addWithDefaultLabel(desc);
0261 }
0262
0263 typedef JetPFConstituentVarProducer<pat::Jet> PatJetPFConstituentVarProducer;
0264 DEFINE_FWK_MODULE(PatJetPFConstituentVarProducer);