File indexing completed on 2024-04-06 12:20:10
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025 #include <memory>
0026
0027 #include <algorithm>
0028 #include <cmath>
0029
0030
0031 #include "FWCore/Framework/interface/Frameworkfwd.h"
0032 #include "FWCore/Framework/interface/stream/EDProducer.h"
0033
0034 #include "FWCore/Framework/interface/Event.h"
0035 #include "FWCore/Framework/interface/MakerMacros.h"
0036 #include "DataFormats/Common/interface/View.h"
0037 #include "DataFormats/Candidate/interface/Candidate.h"
0038 #include "DataFormats/JetReco/interface/CaloJet.h"
0039
0040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0041 #include "FWCore/Utilities/interface/StreamID.h"
0042
0043 class Phase1L1TJetCalibrator : public edm::stream::EDProducer<> {
0044 public:
0045 explicit Phase1L1TJetCalibrator(const edm::ParameterSet&);
0046 ~Phase1L1TJetCalibrator() override;
0047
0048 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0049
0050 private:
0051 void produce(edm::Event&, const edm::EventSetup&) override;
0052
0053 edm::EDGetTokenT<std::vector<reco::CaloJet>> inputCollectionTag_;
0054 std::vector<double> absEtaBinning_;
0055 size_t nBinsEta_;
0056 std::vector<edm::ParameterSet> calibration_;
0057 std::string outputCollectionName_;
0058
0059 std::vector<std::vector<double>> jetCalibrationFactorsBinnedInEta_;
0060 std::vector<std::vector<double>> _jetCalibrationFactorsPtBins;
0061
0062
0063 };
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076 Phase1L1TJetCalibrator::Phase1L1TJetCalibrator(const edm::ParameterSet& iConfig)
0077 : inputCollectionTag_(edm::EDGetTokenT<std::vector<reco::CaloJet>>(
0078 consumes<std::vector<reco::CaloJet>>(iConfig.getParameter<edm::InputTag>("inputCollectionTag")))),
0079 absEtaBinning_(iConfig.getParameter<std::vector<double>>("absEtaBinning")),
0080 nBinsEta_(absEtaBinning_.size() - 1),
0081 calibration_(iConfig.getParameter<std::vector<edm::ParameterSet>>("calibration")),
0082 outputCollectionName_(iConfig.getParameter<std::string>("outputCollectionName"))
0083
0084 {
0085 for (const auto& pset : calibration_) {
0086 _jetCalibrationFactorsPtBins.emplace_back(pset.getParameter<std::vector<double>>("l1tPtBins"));
0087 jetCalibrationFactorsBinnedInEta_.emplace_back(pset.getParameter<std::vector<double>>("l1tCalibrationFactors"));
0088 }
0089
0090 produces<std::vector<reco::CaloJet>>(outputCollectionName_).setBranchAlias(outputCollectionName_);
0091 }
0092
0093 Phase1L1TJetCalibrator::~Phase1L1TJetCalibrator() {}
0094
0095
0096
0097
0098
0099
0100 void Phase1L1TJetCalibrator::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0101 edm::Handle<std::vector<reco::CaloJet>> inputCollectionHandle;
0102 iEvent.getByToken(inputCollectionTag_, inputCollectionHandle);
0103
0104 auto calibratedCollectionPtr = std::make_unique<std::vector<reco::CaloJet>>();
0105
0106
0107
0108
0109
0110
0111
0112
0113
0114
0115 calibratedCollectionPtr->reserve(inputCollectionHandle->size());
0116
0117 for (const auto& candidate : *inputCollectionHandle) {
0118
0119 float pt = candidate.pt();
0120 float eta = candidate.eta();
0121
0122
0123 auto etaBin = upper_bound(absEtaBinning_.begin(), absEtaBinning_.end(), fabs(eta));
0124 int etaIndex = etaBin - absEtaBinning_.begin() - 1;
0125
0126
0127 const std::vector<double>& l1tPtBins = _jetCalibrationFactorsPtBins[etaIndex];
0128 const std::vector<double>& l1tCalibrationFactors = jetCalibrationFactorsBinnedInEta_[etaIndex];
0129
0130
0131 auto ptBin = upper_bound(l1tPtBins.begin(), l1tPtBins.end(), pt);
0132 int ptBinIndex = ptBin - l1tPtBins.begin() - 1;
0133
0134
0135 const double& l1tCalibrationFactor = l1tCalibrationFactors[ptBinIndex];
0136
0137
0138 reco::Candidate::PolarLorentzVector candidateP4(candidate.polarP4());
0139 candidateP4.SetPt(candidateP4.pt() * l1tCalibrationFactor);
0140 calibratedCollectionPtr->emplace_back(candidate);
0141 calibratedCollectionPtr->back().setP4(candidateP4);
0142
0143 #ifdef DEBUG
0144 if (newCandidate->pt() < 0) {
0145 LogDebug("Phase1L1TJetCalibrator") << "######################" << std::endl;
0146 LogDebug("Phase1L1TJetCalibrator") << "PRE-CALIBRATION " << std::endl;
0147 LogDebug("Phase1L1TJetCalibrator") << "\t Jet properties (pt, eta, phi, pile-up): " << candidate.pt() << "\t"
0148 << candidate.eta() << "\t" LogDebug("Phase1L1TJetCalibrator")
0149 << candidate.phi() << "\t" << candidate.pileup() << std::endl;
0150 LogDebug("Phase1L1TJetCalibrator") << "CALIBRATION " << std::endl;
0151 LogDebug("Phase1L1TJetCalibrator") << "\t Using eta - pt - factor " << *etaBin << " - " << *ptBin << " - "
0152 << l1tCalibrationFactor LogDebug("Phase1L1TJetCalibrator") << std::endl;
0153 LogDebug("Phase1L1TJetCalibrator") << "POST-CALIBRATION " << std::endl;
0154 LogDebug("Phase1L1TJetCalibrator") << "\t Jet properties (pt, eta, phi, pile-up): " << newCandidate->pt() << "\t"
0155 << newCandidate->eta() << "\t" << newCandidate->phi() << "\t"
0156 << newCandidate->pileup() << std::endl;
0157 }
0158 #endif
0159 }
0160
0161
0162 std::sort(calibratedCollectionPtr->begin(),
0163 calibratedCollectionPtr->end(),
0164 [](const reco::CaloJet& jet1, const reco::CaloJet& jet2) { return jet1.pt() > jet2.pt(); });
0165
0166 iEvent.put(std::move(calibratedCollectionPtr), outputCollectionName_);
0167 }
0168
0169 void Phase1L1TJetCalibrator::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0170 edm::ParameterSetDescription desc;
0171 desc.add<edm::InputTag>("inputCollectionTag",
0172 edm::InputTag("l1tPhase1JetProducer", "UncalibratedPhase1L1TJetFromPfCandidates"));
0173 desc.add<std::vector<double>>("absEtaBinning");
0174 std::vector<edm::ParameterSet> vDefaults;
0175 edm::ParameterSetDescription validator;
0176 validator.add<double>("etaMax");
0177 validator.add<double>("etaMin");
0178 validator.add<std::vector<double>>("l1tCalibrationFactors");
0179 validator.add<std::vector<double>>("l1tPtBins");
0180 desc.addVPSet("calibration", validator, vDefaults);
0181 desc.add<std::string>("outputCollectionName", "Phase1L1TJetFromPfCandidates");
0182 descriptions.add("Phase1L1TJetCalibrator", desc);
0183 }
0184
0185
0186 DEFINE_FWK_MODULE(Phase1L1TJetCalibrator);