File indexing completed on 2023-03-17 11:17:53
0001 #include "CommonTools/CandAlgos/interface/ModifyObjectValueBase.h"
0002
0003 #include "FWCore/Utilities/interface/InputTag.h"
0004 #include "FWCore/Utilities/interface/EDGetToken.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007
0008 #include "DataFormats/Common/interface/ValueMap.h"
0009 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0010 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0011 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0012 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0013
0014 #include "RecoEgamma/EgammaTools/interface/EgammaRegressionContainer.h"
0015 #include "RecoEgamma/EgammaTools/interface/EpCombinationTool.h"
0016 #include "RecoEgamma/EgammaTools/interface/EcalClusterLocal.h"
0017 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0018 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0019
0020 #include <vdt/vdtMath.h>
0021
0022 class EGRegressionModifierV3 : public ModifyObjectValueBase {
0023 public:
0024 struct EleRegs {
0025 EleRegs(const edm::ParameterSet& iConfig, edm::ConsumesCollector& cc);
0026 void setEventContent(const edm::EventSetup& iSetup);
0027 EgammaRegressionContainer ecalOnlyMean;
0028 EgammaRegressionContainer ecalOnlySigma;
0029 EpCombinationTool epComb;
0030 };
0031
0032 struct PhoRegs {
0033 PhoRegs(const edm::ParameterSet& iConfig, edm::ConsumesCollector& cc);
0034 void setEventContent(const edm::EventSetup& iSetup);
0035 EgammaRegressionContainer ecalOnlyMean;
0036 EgammaRegressionContainer ecalOnlySigma;
0037 };
0038
0039 EGRegressionModifierV3(const edm::ParameterSet& conf, edm::ConsumesCollector& cc);
0040 ~EGRegressionModifierV3() override;
0041
0042 void setEvent(const edm::Event&) final;
0043 void setEventContent(const edm::EventSetup&) final;
0044
0045 void modifyObject(reco::GsfElectron&) const final;
0046 void modifyObject(reco::Photon&) const final;
0047
0048
0049 void modifyObject(pat::Electron&) const final;
0050 void modifyObject(pat::Photon&) const final;
0051
0052 private:
0053 std::array<float, 32> getRegData(const reco::GsfElectron& ele) const;
0054 std::array<float, 32> getRegData(const reco::Photon& pho) const;
0055 void getSeedCrysCoord(const reco::CaloCluster& clus, int& iEtaOrX, int& iPhiOrY) const;
0056
0057 std::unique_ptr<EleRegs> eleRegs_;
0058 std::unique_ptr<PhoRegs> phoRegs_;
0059
0060 float rhoValue_;
0061 edm::EDGetTokenT<double> rhoToken_;
0062
0063 bool useClosestToCentreSeedCrysDef_;
0064 float maxRawEnergyForLowPtEBSigma_;
0065 float maxRawEnergyForLowPtEESigma_;
0066 edm::ESGetToken<CaloGeometry, CaloGeometryRecord> caloGeomToken_;
0067 edm::ESHandle<CaloGeometry> caloGeomHandle_;
0068 };
0069
0070 DEFINE_EDM_PLUGIN(ModifyObjectValueFactory, EGRegressionModifierV3, "EGRegressionModifierV3");
0071
0072 EGRegressionModifierV3::EGRegressionModifierV3(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0073 : ModifyObjectValueBase(conf),
0074 rhoValue_(0.),
0075 rhoToken_(cc.consumes(conf.getParameter<edm::InputTag>("rhoTag"))),
0076 useClosestToCentreSeedCrysDef_(conf.getParameter<bool>("useClosestToCentreSeedCrysDef")),
0077 maxRawEnergyForLowPtEBSigma_(conf.getParameter<double>("maxRawEnergyForLowPtEBSigma")),
0078 maxRawEnergyForLowPtEESigma_(conf.getParameter<double>("maxRawEnergyForLowPtEESigma")) {
0079 if (conf.exists("eleRegs")) {
0080 eleRegs_ = std::make_unique<EleRegs>(conf.getParameterSet("eleRegs"), cc);
0081 }
0082 if (conf.exists("phoRegs")) {
0083 phoRegs_ = std::make_unique<PhoRegs>(conf.getParameterSet("phoRegs"), cc);
0084 }
0085 if (useClosestToCentreSeedCrysDef_) {
0086 caloGeomToken_ = cc.esConsumes();
0087 }
0088 }
0089
0090 EGRegressionModifierV3::~EGRegressionModifierV3() {}
0091
0092 void EGRegressionModifierV3::setEvent(const edm::Event& evt) { rhoValue_ = evt.get(rhoToken_); }
0093
0094 void EGRegressionModifierV3::setEventContent(const edm::EventSetup& iSetup) {
0095 if (eleRegs_)
0096 eleRegs_->setEventContent(iSetup);
0097 if (phoRegs_)
0098 phoRegs_->setEventContent(iSetup);
0099 if (useClosestToCentreSeedCrysDef_) {
0100 caloGeomHandle_ = iSetup.getHandle(caloGeomToken_);
0101 }
0102 }
0103
0104 void EGRegressionModifierV3::modifyObject(reco::GsfElectron& ele) const {
0105
0106
0107 if (!eleRegs_)
0108 return;
0109
0110 const reco::SuperClusterRef& superClus = ele.superCluster();
0111
0112
0113 if (superClus->seed()->seed().det() == DetId::Forward)
0114 return;
0115
0116
0117 bool rescaleDependentValues = superClus->clusters().isAvailable();
0118
0119
0120
0121 if (ele.fbrem() == reco::GsfElectron::ClassificationVariables::kDefaultValue)
0122 return;
0123
0124 auto regData = getRegData(ele);
0125 const float rawEnergy = superClus->rawEnergy();
0126 const float rawESEnergy = superClus->preshowerEnergy();
0127
0128 const float rawEt = rawEnergy * superClus->position().rho() / superClus->position().r();
0129 const bool isSaturated = ele.nSaturatedXtals() != 0;
0130
0131 const float ecalMean = eleRegs_->ecalOnlyMean(rawEt, ele.isEB(), isSaturated, regData.data());
0132 const float ecalMeanCorr = ecalMean > 0. ? ecalMean : 1.0;
0133
0134
0135
0136
0137 if (ele.isEB() && maxRawEnergyForLowPtEBSigma_ >= 0 && eleRegs_->ecalOnlySigma.useLowEtBin(rawEt, isSaturated)) {
0138 regData[0] = std::min(regData[0], maxRawEnergyForLowPtEBSigma_);
0139 }
0140 if (!ele.isEB() && maxRawEnergyForLowPtEESigma_ >= 0 && eleRegs_->ecalOnlySigma.useLowEtBin(rawEt, isSaturated)) {
0141 regData[0] = std::min(regData[0], maxRawEnergyForLowPtEESigma_);
0142 }
0143 const float ecalSigma = eleRegs_->ecalOnlySigma(rawEt, ele.isEB(), isSaturated, regData.data());
0144
0145 const float corrEnergy = (rawEnergy + rawESEnergy) * ecalMeanCorr;
0146 const float corrEnergyErr = corrEnergy * ecalSigma;
0147
0148 ele.setCorrectedEcalEnergy(corrEnergy, rescaleDependentValues);
0149 ele.setCorrectedEcalEnergyError(corrEnergyErr);
0150
0151 std::pair<float, float> combEnergyAndErr = eleRegs_->epComb.combine(ele);
0152 const math::XYZTLorentzVector newP4 = ele.p4() * combEnergyAndErr.first / ele.p4().t();
0153 ele.correctMomentum(newP4, ele.trackMomentumError(), combEnergyAndErr.second);
0154 }
0155
0156 void EGRegressionModifierV3::modifyObject(pat::Electron& ele) const {
0157 modifyObject(static_cast<reco::GsfElectron&>(ele));
0158 }
0159
0160 void EGRegressionModifierV3::modifyObject(reco::Photon& pho) const {
0161
0162
0163 if (!phoRegs_)
0164 return;
0165
0166 const reco::SuperClusterRef& superClus = pho.superCluster();
0167
0168
0169 if (superClus->seed()->seed().det() == DetId::Forward)
0170 return;
0171
0172
0173 if (!superClus->clusters().isAvailable())
0174 return;
0175
0176 auto regData = getRegData(pho);
0177
0178 const float rawEnergy = superClus->rawEnergy();
0179 const float rawESEnergy = superClus->preshowerEnergy();
0180
0181 const float rawEt = rawEnergy * superClus->position().rho() / superClus->position().r();
0182 const bool isSaturated = pho.nSaturatedXtals();
0183 const float ecalMean = phoRegs_->ecalOnlyMean(rawEt, pho.isEB(), isSaturated, regData.data());
0184 const float ecalMeanCorr = ecalMean > 0. ? ecalMean : 1.0;
0185
0186
0187 if (pho.isEB() && maxRawEnergyForLowPtEBSigma_ >= 0 && phoRegs_->ecalOnlySigma.useLowEtBin(rawEt, isSaturated)) {
0188 regData[0] = std::min(regData[0], maxRawEnergyForLowPtEBSigma_);
0189 }
0190 if (!pho.isEB() && maxRawEnergyForLowPtEESigma_ >= 0 && phoRegs_->ecalOnlySigma.useLowEtBin(rawEt, isSaturated)) {
0191 regData[0] = std::min(regData[0], maxRawEnergyForLowPtEESigma_);
0192 }
0193 const float ecalSigma = phoRegs_->ecalOnlySigma(rawEt, pho.isEB(), isSaturated, regData.data());
0194
0195 const double corrEnergy = (rawEnergy + rawESEnergy) * ecalMeanCorr;
0196 const double corrEnergyErr = corrEnergy * ecalSigma;
0197
0198 pho.setCorrectedEnergy(reco::Photon::P4type::regression2, corrEnergy, corrEnergyErr, true);
0199 }
0200
0201 void EGRegressionModifierV3::modifyObject(pat::Photon& pho) const { modifyObject(static_cast<reco::Photon&>(pho)); }
0202
0203 std::array<float, 32> EGRegressionModifierV3::getRegData(const reco::GsfElectron& ele) const {
0204 std::array<float, 32> data;
0205
0206 const reco::SuperClusterRef& superClus = ele.superCluster();
0207 const edm::Ptr<reco::CaloCluster>& seedClus = superClus->seed();
0208
0209 const bool isEB = ele.isEB();
0210 const double rawEnergy = superClus->rawEnergy();
0211 const double rawESEnergy = superClus->preshowerEnergy();
0212 const int numberOfClusters = superClus->clusters().size();
0213 const auto& ssFull5x5 = ele.full5x5_showerShape();
0214
0215 float e5x5Inverse = ssFull5x5.e5x5 != 0. ? vdt::fast_inv(ssFull5x5.e5x5) : 0.;
0216
0217 data[0] = rawEnergy;
0218 data[1] = superClus->etaWidth();
0219 data[2] = superClus->phiWidth();
0220 data[3] = superClus->seed()->energy() / rawEnergy;
0221 data[4] = ssFull5x5.e5x5 / rawEnergy;
0222 data[5] = ele.hcalOverEcalBc();
0223 data[6] = rhoValue_;
0224 data[7] = seedClus->eta() - superClus->position().Eta();
0225 data[8] = reco::deltaPhi(seedClus->phi(), superClus->position().Phi());
0226 data[9] = ssFull5x5.r9;
0227 data[10] = ssFull5x5.sigmaIetaIeta;
0228 data[11] = ssFull5x5.sigmaIetaIphi;
0229 data[12] = ssFull5x5.sigmaIphiIphi;
0230 data[13] = ssFull5x5.eMax * e5x5Inverse;
0231 data[14] = ssFull5x5.e2nd * e5x5Inverse;
0232 data[15] = ssFull5x5.eTop * e5x5Inverse;
0233 data[16] = ssFull5x5.eBottom * e5x5Inverse;
0234 data[17] = ssFull5x5.eLeft * e5x5Inverse;
0235 data[18] = ssFull5x5.eRight * e5x5Inverse;
0236 data[19] = ssFull5x5.e2x5Max * e5x5Inverse;
0237 data[20] = ssFull5x5.e2x5Left * e5x5Inverse;
0238 data[21] = ssFull5x5.e2x5Right * e5x5Inverse;
0239 data[22] = ssFull5x5.e2x5Top * e5x5Inverse;
0240 data[23] = ssFull5x5.e2x5Bottom * e5x5Inverse;
0241 data[24] = ele.nSaturatedXtals();
0242 data[25] = std::max(0, numberOfClusters);
0243
0244 if (isEB) {
0245 int iEta, iPhi;
0246 getSeedCrysCoord(*seedClus, iEta, iPhi);
0247 int signIEta = iEta > 0 ? +1 : -1;
0248 data[26] = iEta;
0249 data[27] = iPhi;
0250 data[28] = (iEta - signIEta) % 5;
0251 data[29] = (iPhi - 1) % 2;
0252 const int iEtaCorr = iEta - (iEta > 0 ? +1 : -1);
0253 const int iEtaCorr26 = iEta - (iEta > 0 ? +26 : -26);
0254 data[30] = std::abs(iEta) <= 25 ? iEtaCorr % 20 : iEtaCorr26 % 20;
0255 data[31] = (iPhi - 1) % 20;
0256 } else {
0257 int iX, iY;
0258 getSeedCrysCoord(*seedClus, iX, iY);
0259 data[26] = iX;
0260 data[27] = iY;
0261 data[28] = rawESEnergy / rawEnergy;
0262 }
0263
0264 return data;
0265 }
0266
0267 std::array<float, 32> EGRegressionModifierV3::getRegData(const reco::Photon& pho) const {
0268 std::array<float, 32> data;
0269
0270 const reco::SuperClusterRef& superClus = pho.superCluster();
0271 const edm::Ptr<reco::CaloCluster>& seedClus = superClus->seed();
0272
0273 const bool isEB = pho.isEB();
0274 const double rawEnergy = superClus->rawEnergy();
0275 const double rawESEnergy = superClus->preshowerEnergy();
0276 const int numberOfClusters = superClus->clusters().size();
0277 const auto& ssFull5x5 = pho.full5x5_showerShapeVariables();
0278
0279 float e5x5Inverse = ssFull5x5.e5x5 != 0. ? vdt::fast_inv(ssFull5x5.e5x5) : 0.;
0280
0281 data[0] = rawEnergy;
0282 data[1] = superClus->etaWidth();
0283 data[2] = superClus->phiWidth();
0284 data[3] = superClus->seed()->energy() / rawEnergy;
0285 data[4] = ssFull5x5.e5x5 / rawEnergy;
0286
0287
0288 data[5] = pho.hadronicOverEm();
0289 data[6] = rhoValue_;
0290 data[7] = seedClus->eta() - superClus->position().Eta();
0291 data[8] = reco::deltaPhi(seedClus->phi(), superClus->position().Phi());
0292 data[9] = pho.full5x5_r9();
0293 data[10] = ssFull5x5.sigmaIetaIeta;
0294
0295
0296 data[11] = ssFull5x5.sigmaIetaIphi;
0297 data[12] = ssFull5x5.sigmaIphiIphi;
0298 data[13] = ssFull5x5.maxEnergyXtal * e5x5Inverse;
0299 data[14] = ssFull5x5.e2nd * e5x5Inverse;
0300 data[15] = ssFull5x5.eTop * e5x5Inverse;
0301 data[16] = ssFull5x5.eBottom * e5x5Inverse;
0302 data[17] = ssFull5x5.eLeft * e5x5Inverse;
0303 data[18] = ssFull5x5.eRight * e5x5Inverse;
0304 data[19] = ssFull5x5.e2x5Max * e5x5Inverse;
0305 data[20] = ssFull5x5.e2x5Left * e5x5Inverse;
0306 data[21] = ssFull5x5.e2x5Right * e5x5Inverse;
0307 data[22] = ssFull5x5.e2x5Top * e5x5Inverse;
0308 data[23] = ssFull5x5.e2x5Bottom * e5x5Inverse;
0309 data[24] = pho.nSaturatedXtals();
0310 data[25] = std::max(0, numberOfClusters);
0311
0312 if (isEB) {
0313 int iEta, iPhi;
0314 getSeedCrysCoord(*seedClus, iEta, iPhi);
0315 data[26] = iEta;
0316 data[27] = iPhi;
0317 int signIEta = iEta > 0 ? +1 : -1;
0318 data[28] = (iEta - signIEta) % 5;
0319 data[29] = (iPhi - 1) % 2;
0320 const int iEtaCorr = iEta - (iEta > 0 ? +1 : -1);
0321 const int iEtaCorr26 = iEta - (iEta > 0 ? +26 : -26);
0322 data[30] = std::abs(iEta) <= 25 ? iEtaCorr % 20 : iEtaCorr26 % 20;
0323 data[31] = (iPhi - 1) % 20;
0324 } else {
0325 int iX, iY;
0326 getSeedCrysCoord(*seedClus, iX, iY);
0327 data[26] = iX;
0328 data[27] = iY;
0329 data[28] = rawESEnergy / rawEnergy;
0330 }
0331
0332 return data;
0333 }
0334
0335 void EGRegressionModifierV3::getSeedCrysCoord(const reco::CaloCluster& clus, int& iEtaOrX, int& iPhiOrY) const {
0336 iEtaOrX = 0;
0337 iPhiOrY = 0;
0338
0339 const bool isEB = clus.seed().subdetId() == EcalBarrel;
0340
0341 if (useClosestToCentreSeedCrysDef_) {
0342 float dummy;
0343 if (isEB) {
0344 egammaTools::localEcalClusterCoordsEB(clus, *caloGeomHandle_, dummy, dummy, iEtaOrX, iPhiOrY, dummy, dummy);
0345 } else {
0346 egammaTools::localEcalClusterCoordsEE(clus, *caloGeomHandle_, dummy, dummy, iEtaOrX, iPhiOrY, dummy, dummy);
0347 }
0348 } else {
0349 if (isEB) {
0350 const EBDetId ebId(clus.seed());
0351 iEtaOrX = ebId.ieta();
0352 iPhiOrY = ebId.iphi();
0353 } else {
0354 const EEDetId eeId(clus.seed());
0355 iEtaOrX = eeId.ix();
0356 iPhiOrY = eeId.iy();
0357 }
0358 }
0359 }
0360
0361 EGRegressionModifierV3::EleRegs::EleRegs(const edm::ParameterSet& iConfig, edm::ConsumesCollector& cc)
0362 : ecalOnlyMean(iConfig.getParameterSet("ecalOnlyMean"), cc),
0363 ecalOnlySigma(iConfig.getParameterSet("ecalOnlySigma"), cc),
0364 epComb(iConfig.getParameterSet("epComb"), std::move(cc)) {}
0365
0366 void EGRegressionModifierV3::EleRegs::setEventContent(const edm::EventSetup& iSetup) {
0367 ecalOnlyMean.setEventContent(iSetup);
0368 ecalOnlySigma.setEventContent(iSetup);
0369 epComb.setEventContent(iSetup);
0370 }
0371
0372 EGRegressionModifierV3::PhoRegs::PhoRegs(const edm::ParameterSet& iConfig, edm::ConsumesCollector& cc)
0373 : ecalOnlyMean(iConfig.getParameterSet("ecalOnlyMean"), cc),
0374 ecalOnlySigma(iConfig.getParameterSet("ecalOnlySigma"), cc) {}
0375
0376 void EGRegressionModifierV3::PhoRegs::setEventContent(const edm::EventSetup& iSetup) {
0377 ecalOnlyMean.setEventContent(iSetup);
0378 ecalOnlySigma.setEventContent(iSetup);
0379 }