File indexing completed on 2024-04-06 12:10:17
0001 #include "EgammaAnalysis/ElectronTools/plugins/RegressionEnergyPatElectronProducer.h"
0002
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/Utilities/interface/EDMException.h"
0011 #include "Geometry/Records/interface/CaloTopologyRecord.h"
0012 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0013 #include "EgammaAnalysis/ElectronTools/interface/SuperClusterHelper.h"
0014
0015 #include <iostream>
0016
0017 RegressionEnergyPatElectronProducer::RegressionEnergyPatElectronProducer(const edm::ParameterSet& cfg) {
0018 inputGsfElectronsToken_ =
0019 mayConsume<reco::GsfElectronCollection>(cfg.getParameter<edm::InputTag>("inputElectronsTag"));
0020 inputPatElectronsToken_ = mayConsume<pat::ElectronCollection>(cfg.getParameter<edm::InputTag>("inputElectronsTag"));
0021 inputCollectionType_ = cfg.getParameter<uint32_t>("inputCollectionType");
0022 rhoInputToken_ = consumes<double>(cfg.getParameter<edm::InputTag>("rhoCollection"));
0023 verticesInputToken_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("vertexCollection"));
0024 energyRegressionType_ = cfg.getParameter<uint32_t>("energyRegressionType");
0025 regressionInputFile_ = cfg.getParameter<std::string>("regressionInputFile");
0026 recHitCollectionEBToken_ = mayConsume<EcalRecHitCollection>(cfg.getParameter<edm::InputTag>("recHitCollectionEB"));
0027 recHitCollectionEEToken_ = mayConsume<EcalRecHitCollection>(cfg.getParameter<edm::InputTag>("recHitCollectionEE"));
0028 ecalTopoToken_ = esConsumes();
0029 caloGeomToken_ = esConsumes();
0030 nameEnergyReg_ = cfg.getParameter<std::string>("nameEnergyReg");
0031 nameEnergyErrorReg_ = cfg.getParameter<std::string>("nameEnergyErrorReg");
0032 debug_ = cfg.getUntrackedParameter<bool>("debug");
0033 useReducedRecHits_ = cfg.getParameter<bool>("useRecHitCollections");
0034 produceValueMaps_ = cfg.getParameter<bool>("produceValueMaps");
0035
0036
0037 if (inputCollectionType_ == 0 && !useReducedRecHits_) {
0038 throw cms::Exception("InconsistentParameters")
0039 << " *** Inconsistent configuration : if you read GsfElectrons, you should set useRecHitCollections to true "
0040 "and provide the correcte values to recHitCollectionEB and recHitCollectionEE (most probably "
0041 "reducedEcalRecHitsEB and reducedEcalRecHitsEE )"
0042 << std::endl;
0043 }
0044
0045 if (inputCollectionType_ == 0 && !produceValueMaps_) {
0046 edm::LogPrint("RegressionEnergyPatElectronProducer")
0047 << " You are running on GsfElectrons and the producer is not configured to produce ValueMaps with the results. "
0048 "In that case, it does not nothing !! ";
0049 }
0050
0051 if (inputCollectionType_ == 0) {
0052
0053 } else if (inputCollectionType_ == 1) {
0054 produces<pat::ElectronCollection>();
0055 } else {
0056 throw cms::Exception("InconsistentParameters")
0057 << " inputCollectionType should be either 0 (GsfElectrons) or 1 (pat::Electrons) ";
0058 }
0059
0060
0061 ElectronEnergyRegressionEvaluate::ElectronEnergyRegressionType type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
0062 if (energyRegressionType_ == 1)
0063 type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
0064 else if (energyRegressionType_ == 2)
0065 type = ElectronEnergyRegressionEvaluate::kWithSubCluVar;
0066 else if (energyRegressionType_ == 3)
0067 type = ElectronEnergyRegressionEvaluate::kWithTrkVar;
0068
0069
0070 regressionEvaluator_ = new ElectronEnergyRegressionEvaluate();
0071 regressionEvaluator_->initialize(regressionInputFile_, type);
0072
0073 if (produceValueMaps_) {
0074 produces<edm::ValueMap<double> >(nameEnergyReg_);
0075 produces<edm::ValueMap<double> >(nameEnergyErrorReg_);
0076 }
0077
0078
0079
0080
0081
0082 geomInitialized_ = false;
0083
0084 edm::LogPrint("RegressionEnergyPatElectronProducer") << " Finished initialization ";
0085 }
0086
0087 RegressionEnergyPatElectronProducer::~RegressionEnergyPatElectronProducer() { delete regressionEvaluator_; }
0088
0089 void RegressionEnergyPatElectronProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0090 assert(regressionEvaluator_->isInitialized());
0091 if (!geomInitialized_) {
0092 ecalTopology_ = &setup.getData(ecalTopoToken_);
0093 caloGeometry_ = &setup.getData(caloGeomToken_);
0094 geomInitialized_ = true;
0095 }
0096
0097
0098
0099
0100 edm::Handle<reco::VertexCollection> hVertexProduct;
0101 event.getByToken(verticesInputToken_, hVertexProduct);
0102 const reco::VertexCollection inVertices = *(hVertexProduct.product());
0103
0104
0105 Int_t nvertices = 0;
0106 for (reco::VertexCollection::const_iterator inV = inVertices.begin(); inV != inVertices.end(); ++inV) {
0107
0108 if (inV->ndof() >= 4 && inV->position().Rho() <= 2.0 && fabs(inV->z()) <= 24.0) {
0109 nvertices++;
0110 }
0111 }
0112
0113
0114
0115
0116 double rho = 0;
0117 edm::Handle<double> hRhoKt6PFJets;
0118 event.getByToken(rhoInputToken_, hRhoKt6PFJets);
0119 rho = (*hRhoKt6PFJets);
0120
0121
0122
0123
0124
0125 edm::Handle<EcalRecHitCollection> pEBRecHits;
0126 edm::Handle<EcalRecHitCollection> pEERecHits;
0127 if (useReducedRecHits_) {
0128 event.getByToken(recHitCollectionEBToken_, pEBRecHits);
0129 event.getByToken(recHitCollectionEEToken_, pEERecHits);
0130 }
0131
0132 edm::Handle<reco::GsfElectronCollection> gsfCollectionH;
0133 edm::Handle<pat::ElectronCollection> patCollectionH;
0134 if (inputCollectionType_ == 0) {
0135 event.getByToken(inputGsfElectronsToken_, gsfCollectionH);
0136 nElectrons_ = gsfCollectionH->size();
0137 }
0138 if (inputCollectionType_ == 1) {
0139 event.getByToken(inputPatElectronsToken_, patCollectionH);
0140 nElectrons_ = patCollectionH->size();
0141 }
0142
0143
0144 std::unique_ptr<pat::ElectronCollection> patElectrons(new pat::ElectronCollection);
0145
0146
0147 std::unique_ptr<edm::ValueMap<double> > regrEnergyMap(new edm::ValueMap<double>());
0148 edm::ValueMap<double>::Filler energyFiller(*regrEnergyMap);
0149
0150 std::unique_ptr<edm::ValueMap<double> > regrEnergyErrorMap(new edm::ValueMap<double>());
0151 edm::ValueMap<double>::Filler energyErrorFiller(*regrEnergyErrorMap);
0152
0153
0154 std::vector<double> energyValues;
0155 std::vector<double> energyErrorValues;
0156 energyValues.reserve(nElectrons_);
0157 energyErrorValues.reserve(nElectrons_);
0158
0159 for (unsigned iele = 0; iele < nElectrons_; ++iele) {
0160 const reco::GsfElectron* ele = (inputCollectionType_ == 0) ? &(*gsfCollectionH)[iele] : &(*patCollectionH)[iele];
0161 if (debug_) {
0162 edm::LogPrint("RegressionEnergyPatElectronProducer")
0163 << "***********************************************************************\n";
0164 edm::LogPrint("RegressionEnergyPatElectronProducer")
0165 << "Run Lumi Event: " << event.id().run() << " " << event.luminosityBlock() << " " << event.id().event()
0166 << "\n";
0167 edm::LogPrint("RegressionEnergyPatElectronProducer")
0168 << "Pat Electron : " << ele->pt() << " " << ele->eta() << " " << ele->phi() << "\n";
0169 }
0170
0171 pat::Electron* myPatElectron = (inputCollectionType_ == 0) ? nullptr : new pat::Electron((*patCollectionH)[iele]);
0172
0173 const EcalRecHitCollection* recHits = nullptr;
0174 if (useReducedRecHits_) {
0175 if (ele->isEB()) {
0176 recHits = pEBRecHits.product();
0177 } else
0178 recHits = pEERecHits.product();
0179 } else {
0180 recHits = (*patCollectionH)[iele].recHits();
0181 }
0182
0183 SuperClusterHelper* mySCHelper = nullptr;
0184 if (inputCollectionType_ == 0) {
0185 mySCHelper = new SuperClusterHelper(&(*ele), recHits, ecalTopology_, caloGeometry_);
0186 } else if (inputCollectionType_ == 1) {
0187 mySCHelper = new SuperClusterHelper(&(*patCollectionH)[iele], recHits, ecalTopology_, caloGeometry_);
0188 }
0189
0190
0191 Double_t FinalMomentum = 0;
0192 Double_t FinalMomentumError = 0;
0193 Double_t RegressionMomentum = 0;
0194 Double_t RegressionMomentumError = 0;
0195
0196 if (energyRegressionType_ == 1) {
0197 RegressionMomentum = regressionEvaluator_->regressionValueNoTrkVar(mySCHelper->rawEnergy(),
0198 mySCHelper->eta(),
0199 mySCHelper->phi(),
0200 mySCHelper->r9(),
0201 mySCHelper->etaWidth(),
0202 mySCHelper->phiWidth(),
0203 mySCHelper->clustersSize(),
0204 mySCHelper->hadronicOverEm(),
0205 rho,
0206 nvertices,
0207 mySCHelper->seedEta(),
0208 mySCHelper->seedPhi(),
0209 mySCHelper->seedEnergy(),
0210 mySCHelper->e3x3(),
0211 mySCHelper->e5x5(),
0212 mySCHelper->sigmaIetaIeta(),
0213 mySCHelper->spp(),
0214 mySCHelper->sep(),
0215 mySCHelper->eMax(),
0216 mySCHelper->e2nd(),
0217 mySCHelper->eTop(),
0218 mySCHelper->eBottom(),
0219 mySCHelper->eLeft(),
0220 mySCHelper->eRight(),
0221 mySCHelper->e2x5Max(),
0222 mySCHelper->e2x5Top(),
0223 mySCHelper->e2x5Bottom(),
0224 mySCHelper->e2x5Left(),
0225 mySCHelper->e2x5Right(),
0226 mySCHelper->ietaSeed(),
0227 mySCHelper->iphiSeed(),
0228 mySCHelper->etaCrySeed(),
0229 mySCHelper->phiCrySeed(),
0230 mySCHelper->preshowerEnergyOverRaw(),
0231 debug_);
0232 RegressionMomentumError =
0233 regressionEvaluator_->regressionUncertaintyNoTrkVar(mySCHelper->rawEnergy(),
0234 mySCHelper->eta(),
0235 mySCHelper->phi(),
0236 mySCHelper->r9(),
0237 mySCHelper->etaWidth(),
0238 mySCHelper->phiWidth(),
0239 mySCHelper->clustersSize(),
0240 mySCHelper->hadronicOverEm(),
0241 rho,
0242 nvertices,
0243 mySCHelper->seedEta(),
0244 mySCHelper->seedPhi(),
0245 mySCHelper->seedEnergy(),
0246 mySCHelper->e3x3(),
0247 mySCHelper->e5x5(),
0248 mySCHelper->sigmaIetaIeta(),
0249 mySCHelper->spp(),
0250 mySCHelper->sep(),
0251 mySCHelper->eMax(),
0252 mySCHelper->e2nd(),
0253 mySCHelper->eTop(),
0254 mySCHelper->eBottom(),
0255 mySCHelper->eLeft(),
0256 mySCHelper->eRight(),
0257 mySCHelper->e2x5Max(),
0258 mySCHelper->e2x5Top(),
0259 mySCHelper->e2x5Bottom(),
0260 mySCHelper->e2x5Left(),
0261 mySCHelper->e2x5Right(),
0262 mySCHelper->ietaSeed(),
0263 mySCHelper->iphiSeed(),
0264 mySCHelper->etaCrySeed(),
0265 mySCHelper->phiCrySeed(),
0266 mySCHelper->preshowerEnergyOverRaw(),
0267 debug_);
0268
0269
0270 if (inputCollectionType_ == 1) {
0271 myPatElectron->setEcalRegressionEnergy(RegressionMomentum, RegressionMomentumError);
0272 }
0273 energyValues.push_back(RegressionMomentum);
0274 energyErrorValues.push_back(RegressionMomentumError);
0275
0276 } else if (energyRegressionType_ == 2) {
0277 RegressionMomentum = regressionEvaluator_->regressionValueWithSubClusters(mySCHelper->rawEnergy(),
0278 mySCHelper->eta(),
0279 mySCHelper->phi(),
0280 mySCHelper->r9(),
0281 mySCHelper->etaWidth(),
0282 mySCHelper->phiWidth(),
0283 mySCHelper->clustersSize(),
0284 mySCHelper->hadronicOverEm(),
0285 rho,
0286 nvertices,
0287 mySCHelper->seedEta(),
0288 mySCHelper->seedPhi(),
0289 mySCHelper->seedEnergy(),
0290 mySCHelper->e3x3(),
0291 mySCHelper->e5x5(),
0292 mySCHelper->sigmaIetaIeta(),
0293 mySCHelper->spp(),
0294 mySCHelper->sep(),
0295 mySCHelper->eMax(),
0296 mySCHelper->e2nd(),
0297 mySCHelper->eTop(),
0298 mySCHelper->eBottom(),
0299 mySCHelper->eLeft(),
0300 mySCHelper->eRight(),
0301 mySCHelper->e2x5Max(),
0302 mySCHelper->e2x5Top(),
0303 mySCHelper->e2x5Bottom(),
0304 mySCHelper->e2x5Left(),
0305 mySCHelper->e2x5Right(),
0306 mySCHelper->ietaSeed(),
0307 mySCHelper->iphiSeed(),
0308 mySCHelper->etaCrySeed(),
0309 mySCHelper->phiCrySeed(),
0310 mySCHelper->preshowerEnergyOverRaw(),
0311 ele->ecalDrivenSeed(),
0312 ele->isEBEtaGap(),
0313 ele->isEBPhiGap(),
0314 ele->isEEDeeGap(),
0315 mySCHelper->eSubClusters(),
0316 mySCHelper->subClusterEnergy(1),
0317 mySCHelper->subClusterEta(1),
0318 mySCHelper->subClusterPhi(1),
0319 mySCHelper->subClusterEmax(1),
0320 mySCHelper->subClusterE3x3(1),
0321 mySCHelper->subClusterEnergy(2),
0322 mySCHelper->subClusterEta(2),
0323 mySCHelper->subClusterPhi(2),
0324 mySCHelper->subClusterEmax(2),
0325 mySCHelper->subClusterE3x3(2),
0326 mySCHelper->subClusterEnergy(3),
0327 mySCHelper->subClusterEta(3),
0328 mySCHelper->subClusterPhi(3),
0329 mySCHelper->subClusterEmax(3),
0330 mySCHelper->subClusterE3x3(3),
0331 mySCHelper->nPreshowerClusters(),
0332 mySCHelper->eESClusters(),
0333 mySCHelper->esClusterEnergy(0),
0334 mySCHelper->esClusterEta(0),
0335 mySCHelper->esClusterPhi(0),
0336 mySCHelper->esClusterEnergy(1),
0337 mySCHelper->esClusterEta(1),
0338 mySCHelper->esClusterPhi(1),
0339 mySCHelper->esClusterEnergy(2),
0340 mySCHelper->esClusterEta(2),
0341 mySCHelper->esClusterPhi(2),
0342 ele->isEB(),
0343 debug_);
0344 RegressionMomentumError =
0345 regressionEvaluator_->regressionUncertaintyWithSubClusters(mySCHelper->rawEnergy(),
0346 mySCHelper->eta(),
0347 mySCHelper->phi(),
0348 mySCHelper->r9(),
0349 mySCHelper->etaWidth(),
0350 mySCHelper->phiWidth(),
0351 mySCHelper->clustersSize(),
0352 mySCHelper->hadronicOverEm(),
0353 rho,
0354 nvertices,
0355 mySCHelper->seedEta(),
0356 mySCHelper->seedPhi(),
0357 mySCHelper->seedEnergy(),
0358 mySCHelper->e3x3(),
0359 mySCHelper->e5x5(),
0360 mySCHelper->sigmaIetaIeta(),
0361 mySCHelper->spp(),
0362 mySCHelper->sep(),
0363 mySCHelper->eMax(),
0364 mySCHelper->e2nd(),
0365 mySCHelper->eTop(),
0366 mySCHelper->eBottom(),
0367 mySCHelper->eLeft(),
0368 mySCHelper->eRight(),
0369 mySCHelper->e2x5Max(),
0370 mySCHelper->e2x5Top(),
0371 mySCHelper->e2x5Bottom(),
0372 mySCHelper->e2x5Left(),
0373 mySCHelper->e2x5Right(),
0374 mySCHelper->ietaSeed(),
0375 mySCHelper->iphiSeed(),
0376 mySCHelper->etaCrySeed(),
0377 mySCHelper->phiCrySeed(),
0378 mySCHelper->preshowerEnergyOverRaw(),
0379 ele->ecalDrivenSeed(),
0380 ele->isEBEtaGap(),
0381 ele->isEBPhiGap(),
0382 ele->isEEDeeGap(),
0383 mySCHelper->eSubClusters(),
0384 mySCHelper->subClusterEnergy(1),
0385 mySCHelper->subClusterEta(1),
0386 mySCHelper->subClusterPhi(1),
0387 mySCHelper->subClusterEmax(1),
0388 mySCHelper->subClusterE3x3(1),
0389 mySCHelper->subClusterEnergy(2),
0390 mySCHelper->subClusterEta(2),
0391 mySCHelper->subClusterPhi(2),
0392 mySCHelper->subClusterEmax(2),
0393 mySCHelper->subClusterE3x3(2),
0394 mySCHelper->subClusterEnergy(3),
0395 mySCHelper->subClusterEta(3),
0396 mySCHelper->subClusterPhi(3),
0397 mySCHelper->subClusterEmax(3),
0398 mySCHelper->subClusterE3x3(3),
0399 mySCHelper->nPreshowerClusters(),
0400 mySCHelper->eESClusters(),
0401 mySCHelper->esClusterEnergy(0),
0402 mySCHelper->esClusterEta(0),
0403 mySCHelper->esClusterPhi(0),
0404 mySCHelper->esClusterEnergy(1),
0405 mySCHelper->esClusterEta(1),
0406 mySCHelper->esClusterPhi(1),
0407 mySCHelper->esClusterEnergy(2),
0408 mySCHelper->esClusterEta(2),
0409 mySCHelper->esClusterPhi(2),
0410 ele->isEB(),
0411 debug_);
0412
0413
0414 if (inputCollectionType_ == 1) {
0415 myPatElectron->setEcalRegressionEnergy(RegressionMomentum, RegressionMomentumError);
0416 }
0417 energyValues.push_back(RegressionMomentum);
0418 energyErrorValues.push_back(RegressionMomentumError);
0419
0420 }
0421
0422 else if (energyRegressionType_ == 3) {
0423 RegressionMomentum = regressionEvaluator_->regressionValueWithTrkVar(ele->p(),
0424 mySCHelper->rawEnergy(),
0425 mySCHelper->eta(),
0426 mySCHelper->phi(),
0427 mySCHelper->etaWidth(),
0428 mySCHelper->phiWidth(),
0429 mySCHelper->clustersSize(),
0430 mySCHelper->hadronicOverEm(),
0431 mySCHelper->r9(),
0432 rho,
0433 nvertices,
0434 mySCHelper->seedEta(),
0435 mySCHelper->seedPhi(),
0436 mySCHelper->seedEnergy(),
0437 mySCHelper->e3x3(),
0438 mySCHelper->e5x5(),
0439 mySCHelper->sigmaIetaIeta(),
0440 mySCHelper->spp(),
0441 mySCHelper->sep(),
0442 mySCHelper->eMax(),
0443 mySCHelper->e2nd(),
0444 mySCHelper->eTop(),
0445 mySCHelper->eBottom(),
0446 mySCHelper->eLeft(),
0447 mySCHelper->eRight(),
0448 mySCHelper->e2x5Max(),
0449 mySCHelper->e2x5Top(),
0450 mySCHelper->e2x5Bottom(),
0451 mySCHelper->e2x5Left(),
0452 mySCHelper->e2x5Right(),
0453 ele->pt(),
0454 ele->trackMomentumAtVtx().R(),
0455 ele->fbrem(),
0456 ele->charge(),
0457 ele->eSuperClusterOverP(),
0458 mySCHelper->ietaSeed(),
0459 mySCHelper->iphiSeed(),
0460 mySCHelper->etaCrySeed(),
0461 mySCHelper->phiCrySeed(),
0462 mySCHelper->preshowerEnergy(),
0463 debug_);
0464 RegressionMomentumError = regressionEvaluator_->regressionUncertaintyWithTrkVar(ele->p(),
0465 mySCHelper->rawEnergy(),
0466 mySCHelper->eta(),
0467 mySCHelper->phi(),
0468 mySCHelper->etaWidth(),
0469 mySCHelper->phiWidth(),
0470 mySCHelper->clustersSize(),
0471 mySCHelper->hadronicOverEm(),
0472 mySCHelper->r9(),
0473 rho,
0474 nvertices,
0475 mySCHelper->seedEta(),
0476 mySCHelper->seedPhi(),
0477 mySCHelper->seedEnergy(),
0478 mySCHelper->e3x3(),
0479 mySCHelper->e5x5(),
0480 mySCHelper->sigmaIetaIeta(),
0481 mySCHelper->spp(),
0482 mySCHelper->sep(),
0483 mySCHelper->eMax(),
0484 mySCHelper->e2nd(),
0485 mySCHelper->eTop(),
0486 mySCHelper->eBottom(),
0487 mySCHelper->eLeft(),
0488 mySCHelper->eRight(),
0489 mySCHelper->e2x5Max(),
0490 mySCHelper->e2x5Top(),
0491 mySCHelper->e2x5Bottom(),
0492 mySCHelper->e2x5Left(),
0493 mySCHelper->e2x5Right(),
0494 ele->pt(),
0495 ele->trackMomentumAtVtx().R(),
0496 ele->fbrem(),
0497 ele->charge(),
0498 ele->eSuperClusterOverP(),
0499 mySCHelper->ietaSeed(),
0500 mySCHelper->iphiSeed(),
0501 mySCHelper->etaCrySeed(),
0502 mySCHelper->phiCrySeed(),
0503 mySCHelper->preshowerEnergy(),
0504 debug_);
0505 FinalMomentum = RegressionMomentum;
0506 FinalMomentumError = RegressionMomentumError;
0507 const math::XYZTLorentzVector& oldMomentum = ele->p4();
0508 math::XYZTLorentzVector newMomentum = math::XYZTLorentzVector(oldMomentum.x() * FinalMomentum / oldMomentum.t(),
0509 oldMomentum.y() * FinalMomentum / oldMomentum.t(),
0510 oldMomentum.z() * FinalMomentum / oldMomentum.t(),
0511 FinalMomentum);
0512
0513 myPatElectron->correctEcalEnergy(RegressionMomentum, RegressionMomentumError);
0514 myPatElectron->correctMomentum(newMomentum, ele->trackMomentumError(), FinalMomentumError);
0515
0516 energyValues.push_back(RegressionMomentum);
0517 energyErrorValues.push_back(RegressionMomentumError);
0518 } else {
0519 edm::LogPrint("RegressionEnergyPatElectronProducer")
0520 << "Error: RegressionType = " << energyRegressionType_ << " is not supported.\n";
0521 }
0522
0523 if (inputCollectionType_ == 1) {
0524 patElectrons->push_back(*myPatElectron);
0525 }
0526 if (myPatElectron)
0527 delete myPatElectron;
0528 if (mySCHelper)
0529 delete mySCHelper;
0530 }
0531
0532
0533 if (inputCollectionType_ == 1) {
0534 event.put(std::move(patElectrons));
0535 }
0536
0537
0538 if (produceValueMaps_) {
0539 if (inputCollectionType_ == 0) {
0540 energyFiller.insert(gsfCollectionH, energyValues.begin(), energyValues.end());
0541 energyErrorFiller.insert(gsfCollectionH, energyErrorValues.begin(), energyErrorValues.end());
0542 } else if (inputCollectionType_ == 1) {
0543 energyFiller.insert(patCollectionH, energyValues.begin(), energyValues.end());
0544 energyErrorFiller.insert(patCollectionH, energyErrorValues.begin(), energyErrorValues.end());
0545 }
0546
0547 energyFiller.fill();
0548 energyErrorFiller.fill();
0549 event.put(std::move(regrEnergyMap), nameEnergyReg_);
0550 event.put(std::move(regrEnergyErrorMap), nameEnergyErrorReg_);
0551 }
0552 }
0553
0554 #include "FWCore/Framework/interface/MakerMacros.h"
0555 #include "FWCore/Framework/interface/ESProducer.h"
0556 #include "FWCore/Framework/interface/ModuleFactory.h"
0557 DEFINE_FWK_MODULE(RegressionEnergyPatElectronProducer);