File indexing completed on 2024-04-06 12:24:05
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "FWCore/Framework/interface/global/EDProducer.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "DataFormats/Common/interface/ValueMap.h"
0017
0018 #include "DataFormats/METReco/interface/CorrMETData.h"
0019 #include "DataFormats/Common/interface/View.h"
0020 #include "DataFormats/Candidate/interface/Candidate.h"
0021
0022 #include <string>
0023 #include <vector>
0024
0025 class ShiftedParticleMETcorrInputProducer : public edm::global::EDProducer<> {
0026 public:
0027 explicit ShiftedParticleMETcorrInputProducer(const edm::ParameterSet&);
0028 ~ShiftedParticleMETcorrInputProducer() override;
0029
0030 private:
0031 typedef edm::View<reco::Candidate> CandidateView;
0032
0033 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0034
0035 const edm::EDGetTokenT<CandidateView> srcOriginalToken_;
0036 const edm::EDGetTokenT<CandidateView> srcShiftedToken_;
0037 edm::EDGetTokenT<edm::ValueMap<float>> weightsToken_;
0038 };
0039
0040 ShiftedParticleMETcorrInputProducer::ShiftedParticleMETcorrInputProducer(const edm::ParameterSet& cfg)
0041 : srcOriginalToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcOriginal"))),
0042 srcShiftedToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcShifted"))) {
0043 edm::InputTag srcWeights = cfg.getParameter<edm::InputTag>("srcWeights");
0044 if (!srcWeights.label().empty())
0045 weightsToken_ = consumes<edm::ValueMap<float>>(srcWeights);
0046
0047 produces<CorrMETData>();
0048 }
0049
0050 ShiftedParticleMETcorrInputProducer::~ShiftedParticleMETcorrInputProducer() {
0051
0052 }
0053
0054 void ShiftedParticleMETcorrInputProducer::produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& es) const {
0055 edm::Handle<CandidateView> originalParticles;
0056 evt.getByToken(srcOriginalToken_, originalParticles);
0057
0058 edm::Handle<CandidateView> shiftedParticles;
0059 evt.getByToken(srcShiftedToken_, shiftedParticles);
0060
0061 edm::Handle<edm::ValueMap<float>> weights;
0062 if (!weightsToken_.isUninitialized())
0063 evt.getByToken(weightsToken_, weights);
0064
0065 auto metCorrection = std::make_unique<CorrMETData>();
0066 if ((!weightsToken_.isUninitialized()) && (originalParticles->size() != shiftedParticles->size()))
0067 throw cms::Exception("InvalidInput")
0068 << "Original collection and shifted collection are of different size in ShiftedParticleMETcorrInputProducer\n";
0069 for (unsigned i = 0; i < originalParticles->size(); ++i) {
0070 float weight = 1.0;
0071 if (!weightsToken_.isUninitialized()) {
0072 edm::Ptr<reco::Candidate> particlePtr = originalParticles->ptrAt(i);
0073 while (!weights->contains(particlePtr.id()) && (particlePtr->numberOfSourceCandidatePtrs() > 0))
0074 particlePtr = particlePtr->sourceCandidatePtr(0);
0075 weight = (*weights)[particlePtr];
0076 }
0077 const reco::Candidate& originalParticle = originalParticles->at(i);
0078 metCorrection->mex += originalParticle.px() * weight;
0079 metCorrection->mey += originalParticle.py() * weight;
0080 metCorrection->sumet -= originalParticle.et() * weight;
0081 }
0082 for (unsigned i = 0; i < shiftedParticles->size(); ++i) {
0083 float weight = 1.0;
0084 if (!weightsToken_.isUninitialized()) {
0085 edm::Ptr<reco::Candidate> particlePtr = originalParticles->ptrAt(i);
0086 while (!weights->contains(particlePtr.id()) && (particlePtr->numberOfSourceCandidatePtrs() > 0))
0087 particlePtr = particlePtr->sourceCandidatePtr(0);
0088 weight = (*weights)[particlePtr];
0089 }
0090 const reco::Candidate& shiftedParticle = shiftedParticles->at(i);
0091 metCorrection->mex -= shiftedParticle.px() * weight;
0092 metCorrection->mey -= shiftedParticle.py() * weight;
0093 metCorrection->sumet += shiftedParticle.et() * weight;
0094 }
0095
0096 evt.put(std::move(metCorrection));
0097 }
0098
0099 #include "FWCore/Framework/interface/MakerMacros.h"
0100
0101 DEFINE_FWK_MODULE(ShiftedParticleMETcorrInputProducer);