Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-22 22:52:36

0001 /** \class ReducedEGProducer
0002  **  
0003  **  Select subset of electrons and photons from input collections and
0004  **  produced consistently relinked output collections including
0005  **  associated SuperClusters, CaloClusters and ecal RecHits
0006  **
0007  **  \author J.Bendavid (CERN)
0008  **  \edited: K. McDermott(Cornell) : refactored code + out of time photons
0009  ***/
0010 
0011 #include "CommonTools/Egamma/interface/ConversionTools.h"
0012 #include "CommonTools/UtilAlgos/interface/StringCutObjectSelector.h"
0013 #include "CommonTools/Utils/interface/StringToEnumValue.h"
0014 #include "CondFormats/EcalObjects/interface/EcalFunctionParameters.h"
0015 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "DataFormats/Common/interface/ValueMap.h"
0018 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
0019 #include "DataFormats/EgammaCandidates/interface/HIPhotonIsolation.h"
0020 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0021 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
0022 #include "DataFormats/EgammaCandidates/interface/PhotonCore.h"
0023 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0024 #include "DataFormats/EgammaReco/interface/BasicClusterShapeAssociation.h"
0025 #include "DataFormats/EgammaReco/interface/ClusterShape.h"
0026 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0027 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0028 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0029 #include "FWCore/Framework/interface/ESHandle.h"
0030 #include "FWCore/Framework/interface/Event.h"
0031 #include "FWCore/Framework/interface/EventSetup.h"
0032 #include "FWCore/Framework/interface/stream/EDProducer.h"
0033 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0035 #include "FWCore/Utilities/interface/Exception.h"
0036 #include "Geometry/Records/interface/CaloTopologyRecord.h"
0037 #include "Geometry/CaloTopology/interface/CaloTopology.h"
0038 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
0039 #include "RecoEgamma/EgammaIsolationAlgos/interface/EGHcalRecHitSelector.h"
0040 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
0041 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgoRcd.h"
0042 
0043 #include <memory>
0044 #include <unordered_set>
0045 #include <vector>
0046 
0047 class ReducedEGProducer : public edm::stream::EDProducer<> {
0048 public:
0049   ReducedEGProducer(const edm::ParameterSet& ps);
0050 
0051   void beginRun(edm::Run const&, const edm::EventSetup&) final;
0052   void produce(edm::Event& evt, const edm::EventSetup& es) final;
0053 
0054 private:
0055   template <class T>
0056   std::vector<edm::EDPutTokenT<T>> vproduces(std::vector<std::string> const& labels) {
0057     std::vector<edm::EDPutTokenT<T>> putTokens{};
0058     putTokens.reserve(labels.size());
0059     for (const auto& label : labels) {
0060       putTokens.push_back(produces<T>(label));
0061     }
0062     return putTokens;
0063   }
0064 
0065   template <typename T, typename U>
0066   void linkCore(const T& core, U& cores, std::map<T, unsigned int>& coreMap);
0067 
0068   void linkSuperCluster(const reco::SuperClusterRef& superCluster,
0069                         std::map<reco::SuperClusterRef, unsigned int>& superClusterMap,
0070                         reco::SuperClusterCollection& superClusters,
0071                         const bool relink,
0072                         std::unordered_set<unsigned int>& superClusterFullRelinkMap);
0073 
0074   void linkConversions(const reco::ConversionRefVector& convrefs,
0075                        reco::ConversionCollection& conversions,
0076                        std::map<reco::ConversionRef, unsigned int>& conversionMap);
0077 
0078   void linkConversionsByTrackRef(const edm::Handle<reco::ConversionCollection>& conversionHandle,
0079                                  const reco::GsfElectron& gsfElectron,
0080                                  reco::ConversionCollection& conversions,
0081                                  std::map<reco::ConversionRef, unsigned int>& conversionMap);
0082 
0083   void linkConversionsByTrackRef(const edm::Handle<reco::ConversionCollection>& conversionHandle,
0084                                  const reco::SuperCluster& superCluster,
0085                                  reco::ConversionCollection& conversions,
0086                                  std::map<reco::ConversionRef, unsigned int>& conversionMap);
0087 
0088   void linkConversion(const reco::ConversionRef& convref,
0089                       reco::ConversionCollection& conversions,
0090                       std::map<reco::ConversionRef, unsigned int>& conversionMap);
0091 
0092   void linkCaloCluster(const reco::CaloClusterPtr& caloCluster,
0093                        reco::CaloClusterCollection& caloClusters,
0094                        std::map<reco::CaloClusterPtr, unsigned int>& caloClusterMap);
0095 
0096   void linkCaloClusters(const reco::SuperCluster& superCluster,
0097                         reco::CaloClusterCollection& ebeeClusters,
0098                         std::map<reco::CaloClusterPtr, unsigned int>& ebeeClusterMap,
0099                         std::unordered_set<DetId>& rechitMap,
0100                         const edm::Handle<EcalRecHitCollection>& barrelHitHandle,
0101                         const edm::Handle<EcalRecHitCollection>& endcapHitHandle,
0102                         CaloTopology const& caloTopology,
0103                         reco::CaloClusterCollection& esClusters,
0104                         std::map<reco::CaloClusterPtr, unsigned int>& esClusterMap);
0105 
0106   void linkHcalHits(const reco::SuperCluster& superClus,
0107                     const HBHERecHitCollection& recHits,
0108                     std::unordered_set<DetId>& hcalDetIds);
0109 
0110   void relinkCaloClusters(reco::SuperCluster& superCluster,
0111                           const std::map<reco::CaloClusterPtr, unsigned int>& ebeeClusterMap,
0112                           const std::map<reco::CaloClusterPtr, unsigned int>& esClusterMap,
0113                           const edm::OrphanHandle<reco::CaloClusterCollection>& outEBEEClusterHandle,
0114                           const edm::OrphanHandle<reco::CaloClusterCollection>& outESClusterHandle);
0115 
0116   template <typename T>
0117   void relinkSuperCluster(T& core,
0118                           const std::map<reco::SuperClusterRef, unsigned int>& superClusterMap,
0119                           const edm::OrphanHandle<reco::SuperClusterCollection>& outSuperClusterHandle);
0120 
0121   void relinkGsfTrack(reco::GsfElectronCore& electroncore,
0122                       const std::map<reco::GsfTrackRef, unsigned int>& gsfTrackMap,
0123                       const edm::OrphanHandle<reco::GsfTrackCollection>& outGsfTrackHandle);
0124 
0125   void relinkConversions(reco::PhotonCore& photonCore,
0126                          const reco::ConversionRefVector& convrefs,
0127                          const std::map<reco::ConversionRef, unsigned int>& conversionMap,
0128                          const edm::OrphanHandle<reco::ConversionCollection>& outConversionHandle);
0129 
0130   void relinkPhotonCore(reco::Photon& photon,
0131                         const std::map<reco::PhotonCoreRef, unsigned int>& photonCoreMap,
0132                         const edm::OrphanHandle<reco::PhotonCoreCollection>& outPhotonCoreHandle);
0133 
0134   void relinkGsfElectronCore(reco::GsfElectron& gsfElectron,
0135                              const std::map<reco::GsfElectronCoreRef, unsigned int>& gsfElectronCoreMap,
0136                              const edm::OrphanHandle<reco::GsfElectronCoreCollection>& outGsfElectronCoreHandle);
0137 
0138   static void calibratePhoton(reco::Photon& photon,
0139                               const reco::PhotonRef& oldPhoRef,
0140                               const edm::ValueMap<float>& energyMap,
0141                               const edm::ValueMap<float>& energyErrMap);
0142 
0143   static void calibrateElectron(reco::GsfElectron& gsfElectron,
0144                                 const reco::GsfElectronRef& oldEleRef,
0145                                 const edm::ValueMap<float>& energyMap,
0146                                 const edm::ValueMap<float>& energyErrMap,
0147                                 const edm::ValueMap<float>& ecalEnergyMap,
0148                                 const edm::ValueMap<float>& ecalEnergyErrMap);
0149 
0150   template <typename T>
0151   void setToken(edm::EDGetTokenT<T>& token, const edm::ParameterSet& config, const std::string& name) {
0152     token = consumes<T>(config.getParameter<edm::InputTag>(name));
0153   }
0154 
0155   //tokens for input collections
0156   const edm::EDGetTokenT<reco::SuperClusterCollection> superclusterT_;
0157   const edm::EDGetTokenT<reco::PhotonCollection> photonT_;
0158   edm::EDGetTokenT<reco::PhotonCollection> ootPhotonT_;
0159   const edm::EDGetTokenT<reco::GsfElectronCollection> gsfElectronT_;
0160   const edm::EDGetTokenT<reco::ConversionCollection> conversionT_;
0161   const edm::EDGetTokenT<reco::ConversionCollection> singleConversionT_;
0162 
0163   const edm::EDGetTokenT<EcalRecHitCollection> barrelEcalHits_;
0164   const edm::EDGetTokenT<EcalRecHitCollection> endcapEcalHits_;
0165   const bool doPreshowerEcalHits_;
0166   const edm::EDGetTokenT<EcalRecHitCollection> preshowerEcalHits_;
0167   const edm::EDGetTokenT<HBHERecHitCollection> hbheHits_;
0168 
0169   const edm::EDGetTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> photonPfCandMapT_;
0170   const edm::EDGetTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> gsfElectronPfCandMapT_;
0171 
0172   std::vector<edm::EDGetTokenT<edm::ValueMap<bool>>> photonIdTs_;
0173   std::vector<edm::EDGetTokenT<edm::ValueMap<float>>> gsfElectronIdTs_;
0174 
0175   std::vector<edm::EDGetTokenT<edm::ValueMap<float>>> photonFloatValueMapTs_;
0176   std::vector<edm::EDGetTokenT<edm::ValueMap<float>>> ootPhotonFloatValueMapTs_;
0177   std::vector<edm::EDGetTokenT<edm::ValueMap<float>>> gsfElectronFloatValueMapTs_;
0178 
0179   const edm::EDGetTokenT<reco::HIPhotonIsolationMap> recoHIPhotonIsolationMapInputToken_;
0180   edm::EDPutTokenT<reco::HIPhotonIsolationMap> recoHIPhotonIsolationMapOutputName_;
0181 
0182   const double scPtMin_;
0183   const double scAbsetaMax_;
0184   const double relinkSuperclusterPtMin_;
0185 
0186   const bool applyPhotonCalibOnData_;
0187   const bool applyPhotonCalibOnMC_;
0188   const bool applyGsfElectronCalibOnData_;
0189   const bool applyGsfElectronCalibOnMC_;
0190   edm::EDGetTokenT<edm::ValueMap<float>> photonCalibEnergyT_;
0191   edm::EDGetTokenT<edm::ValueMap<float>> photonCalibEnergyErrT_;
0192   edm::EDGetTokenT<edm::ValueMap<float>> gsfElectronCalibEnergyT_;
0193   edm::EDGetTokenT<edm::ValueMap<float>> gsfElectronCalibEnergyErrT_;
0194   edm::EDGetTokenT<edm::ValueMap<float>> gsfElectronCalibEcalEnergyT_;
0195   edm::EDGetTokenT<edm::ValueMap<float>> gsfElectronCalibEcalEnergyErrT_;
0196 
0197   edm::ESGetToken<CaloTopology, CaloTopologyRecord> caloTopology_;
0198   //names for output collections
0199   const edm::EDPutTokenT<reco::PhotonCollection> outPhotons_;
0200   const edm::EDPutTokenT<reco::PhotonCoreCollection> outPhotonCores_;
0201   edm::EDPutTokenT<reco::PhotonCollection> outOOTPhotons_;
0202   edm::EDPutTokenT<reco::PhotonCoreCollection> outOOTPhotonCores_;
0203   const edm::EDPutTokenT<reco::GsfElectronCollection> outGsfElectrons_;
0204   const edm::EDPutTokenT<reco::GsfElectronCoreCollection> outGsfElectronCores_;
0205   const edm::EDPutTokenT<reco::GsfTrackCollection> outGsfTracks_;
0206   const edm::EDPutTokenT<reco::ConversionCollection> outConversions_;
0207   const edm::EDPutTokenT<reco::ConversionCollection> outSingleConversions_;
0208   const edm::EDPutTokenT<reco::SuperClusterCollection> outSuperClusters_;
0209   const edm::EDPutTokenT<reco::CaloClusterCollection> outEBEEClusters_;
0210   const edm::EDPutTokenT<reco::CaloClusterCollection> outESClusters_;
0211   edm::EDPutTokenT<reco::SuperClusterCollection> outOOTSuperClusters_;
0212   edm::EDPutTokenT<reco::CaloClusterCollection> outOOTEBEEClusters_;
0213   edm::EDPutTokenT<reco::CaloClusterCollection> outOOTESClusters_;
0214   const edm::EDPutTokenT<EcalRecHitCollection> outEBRecHits_;
0215   const edm::EDPutTokenT<EcalRecHitCollection> outEERecHits_;
0216   edm::EDPutTokenT<EcalRecHitCollection> outESRecHits_;
0217   const edm::EDPutTokenT<HBHERecHitCollection> outHBHERecHits_;
0218   const edm::EDPutTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> outPhotonPfCandMap_;
0219   const edm::EDPutTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> outGsfElectronPfCandMap_;
0220   const std::vector<edm::EDPutTokenT<edm::ValueMap<bool>>> outPhotonIds_;
0221   const std::vector<edm::EDPutTokenT<edm::ValueMap<float>>> outGsfElectronIds_;
0222   const std::vector<edm::EDPutTokenT<edm::ValueMap<float>>> outPhotonFloatValueMaps_;
0223   std::vector<edm::EDPutTokenT<edm::ValueMap<float>>> outOOTPhotonFloatValueMaps_;
0224   const std::vector<edm::EDPutTokenT<edm::ValueMap<float>>> outGsfElectronFloatValueMaps_;
0225 
0226   const StringCutObjectSelector<reco::Photon> keepPhotonSel_;
0227   const StringCutObjectSelector<reco::Photon> slimRelinkPhotonSel_;
0228   const StringCutObjectSelector<reco::Photon> relinkPhotonSel_;
0229   const StringCutObjectSelector<reco::Photon> keepOOTPhotonSel_;
0230   const StringCutObjectSelector<reco::Photon> slimRelinkOOTPhotonSel_;
0231   const StringCutObjectSelector<reco::Photon> relinkOOTPhotonSel_;
0232   const StringCutObjectSelector<reco::GsfElectron> keepGsfElectronSel_;
0233   const StringCutObjectSelector<reco::GsfElectron> slimRelinkGsfElectronSel_;
0234   const StringCutObjectSelector<reco::GsfElectron> relinkGsfElectronSel_;
0235 
0236   EGHcalRecHitSelector hcalHitSel_;
0237 };
0238 
0239 #include "FWCore/Framework/interface/MakerMacros.h"
0240 DEFINE_FWK_MODULE(ReducedEGProducer);
0241 
0242 namespace {
0243 
0244   template <class T>
0245   auto getHandles(edm::Event const& event, std::vector<edm::EDGetTokenT<T>> const& tokens) {
0246     std::vector<edm::Handle<T>> handles(tokens.size());
0247     int index = 0;
0248     for (const auto& token : tokens) {
0249       event.getByToken(token, handles[index++]);
0250     }
0251     return handles;
0252   }
0253 
0254   template <class Handle, class T>
0255   auto emplaceValueMap(Handle const& handle,
0256                        std::vector<T> const& values,
0257                        edm::Event& ev,
0258                        edm::EDPutTokenT<edm::ValueMap<T>> const& putToken) {
0259     using MapType = edm::ValueMap<T>;
0260     MapType oMap{};
0261     {
0262       typename MapType::Filler filler(oMap);
0263       filler.insert(handle, values.begin(), values.end());
0264       filler.fill();
0265     }
0266     ev.emplace(putToken, std::move(oMap));
0267   };
0268 }  // namespace
0269 
0270 ReducedEGProducer::ReducedEGProducer(const edm::ParameterSet& config)
0271     : superclusterT_(consumes(config.getParameter<edm::InputTag>("pflowSuperclusters"))),
0272       photonT_(consumes(config.getParameter<edm::InputTag>("photons"))),
0273       gsfElectronT_(consumes(config.getParameter<edm::InputTag>("gsfElectrons"))),
0274       conversionT_(consumes(config.getParameter<edm::InputTag>("conversions"))),
0275       singleConversionT_(consumes(config.getParameter<edm::InputTag>("singleConversions"))),
0276       barrelEcalHits_(consumes(config.getParameter<edm::InputTag>("barrelEcalHits"))),
0277       endcapEcalHits_(consumes(config.getParameter<edm::InputTag>("endcapEcalHits"))),
0278       doPreshowerEcalHits_(!config.getParameter<edm::InputTag>("preshowerEcalHits").label().empty()),
0279       preshowerEcalHits_(doPreshowerEcalHits_
0280                              ? consumes<EcalRecHitCollection>(config.getParameter<edm::InputTag>("preshowerEcalHits"))
0281                              : edm::EDGetTokenT<EcalRecHitCollection>()),
0282       hbheHits_(consumes<HBHERecHitCollection>(config.getParameter<edm::InputTag>("hbheHits"))),
0283       photonPfCandMapT_(consumes(config.getParameter<edm::InputTag>("photonsPFValMap"))),
0284       gsfElectronPfCandMapT_(consumes(config.getParameter<edm::InputTag>("gsfElectronsPFValMap"))),
0285       recoHIPhotonIsolationMapInputToken_{
0286           !config.getParameter<edm::InputTag>("hiPhotonIsolationMapInput").label().empty()
0287               ? consumes<reco::HIPhotonIsolationMap>(config.getParameter<edm::InputTag>("hiPhotonIsolationMapInput"))
0288               : edm::EDGetTokenT<reco::HIPhotonIsolationMap>{}},
0289       scPtMin_(config.getParameter<double>("keepPfSuperclusterPtMin")),
0290       scAbsetaMax_(config.getParameter<double>("keepPfSuperclusterAbsetaMax")),
0291       relinkSuperclusterPtMin_(config.getParameter<double>("relinkSuperclusterPtMin")),
0292       //calibration flags
0293       applyPhotonCalibOnData_(config.getParameter<bool>("applyPhotonCalibOnData")),
0294       applyPhotonCalibOnMC_(config.getParameter<bool>("applyPhotonCalibOnMC")),
0295       applyGsfElectronCalibOnData_(config.getParameter<bool>("applyGsfElectronCalibOnData")),
0296       applyGsfElectronCalibOnMC_(config.getParameter<bool>("applyGsfElectronCalibOnMC")),
0297       //output collections
0298       outPhotons_{produces<reco::PhotonCollection>("reducedGedPhotons")},
0299       outPhotonCores_{produces<reco::PhotonCoreCollection>("reducedGedPhotonCores")},
0300       outGsfElectrons_{produces<reco::GsfElectronCollection>("reducedGedGsfElectrons")},
0301       outGsfElectronCores_{produces<reco::GsfElectronCoreCollection>("reducedGedGsfElectronCores")},
0302       outGsfTracks_{produces<reco::GsfTrackCollection>("reducedGsfTracks")},
0303       outConversions_{produces<reco::ConversionCollection>("reducedConversions")},
0304       outSingleConversions_{produces<reco::ConversionCollection>("reducedSingleLegConversions")},
0305       outSuperClusters_{produces<reco::SuperClusterCollection>("reducedSuperClusters")},
0306       outEBEEClusters_{produces<reco::CaloClusterCollection>("reducedEBEEClusters")},
0307       outESClusters_{produces<reco::CaloClusterCollection>("reducedESClusters")},
0308       outEBRecHits_{produces<EcalRecHitCollection>("reducedEBRecHits")},
0309       outEERecHits_{produces<EcalRecHitCollection>("reducedEERecHits")},
0310       outHBHERecHits_{produces<HBHERecHitCollection>("reducedHBHEHits")},
0311       outPhotonPfCandMap_{produces<edm::ValueMap<std::vector<reco::PFCandidateRef>>>("reducedPhotonPfCandMap")},
0312       outGsfElectronPfCandMap_{
0313           produces<edm::ValueMap<std::vector<reco::PFCandidateRef>>>("reducedGsfElectronPfCandMap")},
0314       outPhotonIds_{vproduces<edm::ValueMap<bool>>(config.getParameter<std::vector<std::string>>("photonIDOutput"))},
0315       outGsfElectronIds_{
0316           vproduces<edm::ValueMap<float>>(config.getParameter<std::vector<std::string>>("gsfElectronIDOutput"))},
0317       outPhotonFloatValueMaps_{
0318           vproduces<edm::ValueMap<float>>(config.getParameter<std::vector<std::string>>("photonFloatValueMapOutput"))},
0319       outGsfElectronFloatValueMaps_{vproduces<edm::ValueMap<float>>(
0320           config.getParameter<std::vector<std::string>>("gsfElectronFloatValueMapOutput"))},
0321       keepPhotonSel_(config.getParameter<std::string>("keepPhotons")),
0322       slimRelinkPhotonSel_(config.getParameter<std::string>("slimRelinkPhotons")),
0323       relinkPhotonSel_(config.getParameter<std::string>("relinkPhotons")),
0324       keepOOTPhotonSel_(config.getParameter<std::string>("keepOOTPhotons")),
0325       slimRelinkOOTPhotonSel_(config.getParameter<std::string>("slimRelinkOOTPhotons")),
0326       relinkOOTPhotonSel_(config.getParameter<std::string>("relinkOOTPhotons")),
0327       keepGsfElectronSel_(config.getParameter<std::string>("keepGsfElectrons")),
0328       slimRelinkGsfElectronSel_(config.getParameter<std::string>("slimRelinkGsfElectrons")),
0329       relinkGsfElectronSel_(config.getParameter<std::string>("relinkGsfElectrons")),
0330       hcalHitSel_(config.getParameter<edm::ParameterSet>("hcalHitSel"), consumesCollector()) {
0331   const auto& aTag = config.getParameter<edm::InputTag>("ootPhotons");
0332   caloTopology_ = esConsumes();
0333   if (not aTag.label().empty())
0334     ootPhotonT_ = consumes<reco::PhotonCollection>(aTag);
0335 
0336   for (const edm::InputTag& tag : config.getParameter<std::vector<edm::InputTag>>("photonIDSources")) {
0337     photonIdTs_.emplace_back(consumes<edm::ValueMap<bool>>(tag));
0338   }
0339 
0340   for (const edm::InputTag& tag : config.getParameter<std::vector<edm::InputTag>>("gsfElectronIDSources")) {
0341     gsfElectronIdTs_.emplace_back(consumes<edm::ValueMap<float>>(tag));
0342   }
0343 
0344   for (const edm::InputTag& tag : config.getParameter<std::vector<edm::InputTag>>("photonFloatValueMapSources")) {
0345     photonFloatValueMapTs_.emplace_back(consumes<edm::ValueMap<float>>(tag));
0346   }
0347 
0348   for (const edm::InputTag& tag : config.getParameter<std::vector<edm::InputTag>>("ootPhotonFloatValueMapSources")) {
0349     ootPhotonFloatValueMapTs_.emplace_back(consumes<edm::ValueMap<float>>(tag));
0350   }
0351 
0352   for (const edm::InputTag& tag : config.getParameter<std::vector<edm::InputTag>>("gsfElectronFloatValueMapSources")) {
0353     gsfElectronFloatValueMapTs_.emplace_back(consumes<edm::ValueMap<float>>(tag));
0354   }
0355 
0356   if (applyPhotonCalibOnData_ || applyPhotonCalibOnMC_) {
0357     setToken(photonCalibEnergyT_, config, "photonCalibEnergySource");
0358     setToken(photonCalibEnergyErrT_, config, "photonCalibEnergyErrSource");
0359   }
0360   if (applyGsfElectronCalibOnData_ || applyGsfElectronCalibOnMC_) {
0361     setToken(gsfElectronCalibEnergyT_, config, "gsfElectronCalibEnergySource");
0362     setToken(gsfElectronCalibEnergyErrT_, config, "gsfElectronCalibEnergyErrSource");
0363     setToken(gsfElectronCalibEcalEnergyT_, config, "gsfElectronCalibEcalEnergySource");
0364     setToken(gsfElectronCalibEcalEnergyErrT_, config, "gsfElectronCalibEcalEnergyErrSource");
0365   }
0366 
0367   if (!ootPhotonT_.isUninitialized()) {
0368     outOOTPhotons_ = produces<reco::PhotonCollection>("reducedOOTPhotons");
0369     outOOTPhotonCores_ = produces<reco::PhotonCoreCollection>("reducedOOTPhotonCores");
0370     outOOTSuperClusters_ = produces<reco::SuperClusterCollection>("reducedOOTSuperClusters");
0371     outOOTEBEEClusters_ = produces<reco::CaloClusterCollection>("reducedOOTEBEEClusters");
0372     outOOTESClusters_ = produces<reco::CaloClusterCollection>("reducedOOTESClusters");
0373   }
0374   if (doPreshowerEcalHits_) {
0375     outESRecHits_ = produces<EcalRecHitCollection>("reducedESRecHits");
0376   }
0377   if (!ootPhotonT_.isUninitialized()) {
0378     outOOTPhotonFloatValueMaps_ =
0379         vproduces<edm::ValueMap<float>>(config.getParameter<std::vector<std::string>>("ootPhotonFloatValueMapOutput"));
0380   }
0381   if (!recoHIPhotonIsolationMapInputToken_.isUninitialized()) {
0382     recoHIPhotonIsolationMapOutputName_ =
0383         produces<reco::HIPhotonIsolationMap>(config.getParameter<std::string>("hiPhotonIsolationMapOutput"));
0384   }
0385 }
0386 
0387 void ReducedEGProducer::beginRun(edm::Run const& run, const edm::EventSetup& iSetup) { hcalHitSel_.setup(iSetup); }
0388 
0389 void ReducedEGProducer::produce(edm::Event& event, const edm::EventSetup& eventSetup) {
0390   //get input collections
0391 
0392   auto scHandle = event.getHandle(superclusterT_);
0393   auto photonHandle = event.getHandle(photonT_);
0394 
0395   auto ootPhotonHandle =
0396       !ootPhotonT_.isUninitialized() ? event.getHandle(ootPhotonT_) : edm::Handle<reco::PhotonCollection>{};
0397 
0398   auto gsfElectronHandle = event.getHandle(gsfElectronT_);
0399   auto conversionHandle = event.getHandle(conversionT_);
0400   auto singleConversionHandle = event.getHandle(singleConversionT_);
0401   auto barrelHitHandle = event.getHandle(barrelEcalHits_);
0402   auto endcapHitHandle = event.getHandle(endcapEcalHits_);
0403 
0404   auto preshowerHitHandle =
0405       doPreshowerEcalHits_ ? event.getHandle(preshowerEcalHits_) : edm::Handle<EcalRecHitCollection>{};
0406 
0407   auto hbheHitHandle = event.getHandle(hbheHits_);
0408   auto photonPfCandMapHandle = event.getHandle(photonPfCandMapT_);
0409   auto gsfElectronPfCandMapHandle = event.getHandle(gsfElectronPfCandMapT_);
0410 
0411   auto photonIdHandles = getHandles(event, photonIdTs_);
0412   auto gsfElectronIdHandles = getHandles(event, gsfElectronIdTs_);
0413   auto photonFloatValueMapHandles = getHandles(event, photonFloatValueMapTs_);
0414 
0415   auto ootPhotonFloatValueMapHandles = !ootPhotonT_.isUninitialized()
0416                                            ? getHandles(event, ootPhotonFloatValueMapTs_)
0417                                            : std::vector<edm::Handle<edm::ValueMap<float>>>{};
0418 
0419   auto gsfElectronFloatValueMapHandles = getHandles(event, gsfElectronFloatValueMapTs_);
0420 
0421   edm::Handle<edm::ValueMap<float>> gsfElectronCalibEnergyHandle;
0422   edm::Handle<edm::ValueMap<float>> gsfElectronCalibEnergyErrHandle;
0423   edm::Handle<edm::ValueMap<float>> gsfElectronCalibEcalEnergyHandle;
0424   edm::Handle<edm::ValueMap<float>> gsfElectronCalibEcalEnergyErrHandle;
0425   if (applyGsfElectronCalibOnData_ || applyGsfElectronCalibOnMC_) {
0426     event.getByToken(gsfElectronCalibEnergyT_, gsfElectronCalibEnergyHandle);
0427     event.getByToken(gsfElectronCalibEnergyErrT_, gsfElectronCalibEnergyErrHandle);
0428     event.getByToken(gsfElectronCalibEcalEnergyT_, gsfElectronCalibEcalEnergyHandle);
0429     event.getByToken(gsfElectronCalibEcalEnergyErrT_, gsfElectronCalibEcalEnergyErrHandle);
0430   }
0431   edm::Handle<edm::ValueMap<float>> photonCalibEnergyHandle;
0432   edm::Handle<edm::ValueMap<float>> photonCalibEnergyErrHandle;
0433   if (applyPhotonCalibOnData_ || applyPhotonCalibOnMC_) {
0434     event.getByToken(photonCalibEnergyT_, photonCalibEnergyHandle);
0435     event.getByToken(photonCalibEnergyErrT_, photonCalibEnergyErrHandle);
0436   }
0437 
0438   auto const& caloTopology = eventSetup.getData(caloTopology_);
0439 
0440   //initialize output collections
0441   reco::PhotonCollection photons;
0442   reco::PhotonCoreCollection photonCores;
0443   reco::PhotonCollection ootPhotons;
0444   reco::PhotonCoreCollection ootPhotonCores;
0445   reco::GsfElectronCollection gsfElectrons;
0446   reco::GsfElectronCoreCollection gsfElectronCores;
0447   reco::GsfTrackCollection gsfTracks;
0448   reco::ConversionCollection conversions;
0449   reco::ConversionCollection singleConversions;
0450   reco::SuperClusterCollection superClusters;
0451   reco::CaloClusterCollection ebeeClusters;
0452   reco::CaloClusterCollection esClusters;
0453   reco::SuperClusterCollection ootSuperClusters;
0454   reco::CaloClusterCollection ootEbeeClusters;
0455   reco::CaloClusterCollection ootEsClusters;
0456   EcalRecHitCollection ebRecHits;
0457   EcalRecHitCollection eeRecHits;
0458   EcalRecHitCollection esRecHits;
0459   HBHERecHitCollection hbheRecHits;
0460   edm::ValueMap<std::vector<reco::PFCandidateRef>> photonPfCandMap;
0461   edm::ValueMap<std::vector<reco::PFCandidateRef>> gsfElectronPfCandMap;
0462 
0463   //maps to collection indices of output objects
0464   std::map<reco::PhotonCoreRef, unsigned int> photonCoreMap;
0465   std::map<reco::PhotonCoreRef, unsigned int> ootPhotonCoreMap;
0466   std::map<reco::GsfElectronCoreRef, unsigned int> gsfElectronCoreMap;
0467   std::map<reco::GsfTrackRef, unsigned int> gsfTrackMap;
0468   std::map<reco::ConversionRef, unsigned int> conversionMap;
0469   std::map<reco::ConversionRef, unsigned int> singleConversionMap;
0470   std::map<reco::SuperClusterRef, unsigned int> superClusterMap;
0471   std::map<reco::CaloClusterPtr, unsigned int> ebeeClusterMap;
0472   std::map<reco::CaloClusterPtr, unsigned int> esClusterMap;
0473   std::map<reco::SuperClusterRef, unsigned int> ootSuperClusterMap;
0474   std::map<reco::CaloClusterPtr, unsigned int> ootEbeeClusterMap;
0475   std::map<reco::CaloClusterPtr, unsigned int> ootEsClusterMap;
0476   std::unordered_set<DetId> rechitMap;
0477   std::unordered_set<DetId> hcalRechitMap;
0478 
0479   std::unordered_set<unsigned int> superClusterFullRelinkMap;
0480   std::unordered_set<unsigned int> ootSuperClusterFullRelinkMap;
0481 
0482   //vectors for pfcandidate valuemaps
0483   std::vector<std::vector<reco::PFCandidateRef>> pfCandIsoPairVecPho;
0484   std::vector<std::vector<reco::PFCandidateRef>> pfCandIsoPairVecEle;
0485 
0486   //vectors for id valuemaps
0487   std::vector<std::vector<bool>> photonIdVals(photonIdHandles.size());
0488   std::vector<std::vector<float>> gsfElectronIdVals(gsfElectronIdHandles.size());
0489   std::vector<std::vector<float>> photonFloatValueMapVals(photonFloatValueMapHandles.size());
0490   std::vector<std::vector<float>> ootPhotonFloatValueMapVals(ootPhotonFloatValueMapHandles.size());
0491   std::vector<std::vector<float>> gsfElectronFloatValueMapVals(gsfElectronFloatValueMapHandles.size());
0492 
0493   // HI photon iso value maps
0494   reco::HIPhotonIsolationMap const* recoHIPhotonIsolationMapInputValueMap =
0495       !recoHIPhotonIsolationMapInputToken_.isUninitialized() ? &event.get(recoHIPhotonIsolationMapInputToken_)
0496                                                              : nullptr;
0497   std::vector<reco::HIPhotonIsolation> recoHIPhotonIsolationMapInputVals;
0498 
0499   //loop over photons and fill maps
0500   int index = -1;
0501   for (const auto& photon : *photonHandle) {
0502     index++;
0503 
0504     reco::PhotonRef photonref(photonHandle, index);
0505     photons.push_back(photon);
0506     auto& newPhoton = photons.back();
0507 
0508     if ((applyPhotonCalibOnData_ && event.isRealData()) || (applyPhotonCalibOnMC_ && !event.isRealData())) {
0509       calibratePhoton(newPhoton, photonref, *photonCalibEnergyHandle, *photonCalibEnergyErrHandle);
0510     }
0511 
0512     //we do this after calibration
0513     bool keep = keepPhotonSel_(newPhoton);
0514     if (!keep) {
0515       photons.pop_back();
0516       continue;
0517     }
0518 
0519     //fill pf candidate value map vector
0520     pfCandIsoPairVecPho.push_back((*photonPfCandMapHandle)[photonref]);
0521 
0522     //fill photon id valuemap vectors
0523     int subindex = 0;
0524     for (const auto& photonIdHandle : photonIdHandles) {
0525       photonIdVals[subindex++].push_back((*photonIdHandle)[photonref]);
0526     }
0527 
0528     subindex = 0;
0529     for (const auto& photonFloatValueMapHandle : photonFloatValueMapHandles) {
0530       photonFloatValueMapVals[subindex++].push_back((*photonFloatValueMapHandle)[photonref]);
0531     }
0532 
0533     // HI photon isolation
0534     if (!recoHIPhotonIsolationMapInputToken_.isUninitialized()) {
0535       recoHIPhotonIsolationMapInputVals.push_back((*recoHIPhotonIsolationMapInputValueMap)[photonref]);
0536     }
0537 
0538     //link photon core
0539     const reco::PhotonCoreRef& photonCore = photon.photonCore();
0540     linkCore(photonCore, photonCores, photonCoreMap);
0541 
0542     bool slimRelink = slimRelinkPhotonSel_(newPhoton);
0543     //no supercluster relinking unless slimRelink selection is satisfied
0544     if (!slimRelink)
0545       continue;
0546 
0547     bool relink = relinkPhotonSel_(newPhoton);
0548 
0549     //link supercluster
0550     const reco::SuperClusterRef& superCluster = photon.superCluster();
0551     linkSuperCluster(superCluster, superClusterMap, superClusters, relink, superClusterFullRelinkMap);
0552 
0553     //conversions only for full relinking
0554     if (!relink)
0555       continue;
0556 
0557     const reco::ConversionRefVector& convrefs = photon.conversions();
0558     linkConversions(convrefs, conversions, conversionMap);
0559 
0560     //explicitly references conversions
0561     const reco::ConversionRefVector& singleconvrefs = photon.conversionsOneLeg();
0562     linkConversions(singleconvrefs, singleConversions, singleConversionMap);
0563 
0564     //hcal hits
0565     linkHcalHits(*photon.superCluster(), *hbheHitHandle, hcalRechitMap);
0566   }
0567 
0568   //loop over oot photons and fill maps
0569   //special note1: since not PFCand --> no PF isolation, IDs (but we do have FloatValueMap!)
0570   //special note2: conversion sequence not run over bcs from oot phos, so skip relinking of oot phos
0571   //special note3: clusters and superclusters in own collections!
0572   if (!ootPhotonT_.isUninitialized()) {
0573     index = -1;
0574     for (const auto& ootPhoton : *ootPhotonHandle) {
0575       index++;
0576 
0577       bool keep = keepOOTPhotonSel_(ootPhoton);
0578       if (!keep)
0579         continue;
0580 
0581       reco::PhotonRef ootPhotonref(ootPhotonHandle, index);
0582 
0583       ootPhotons.push_back(ootPhoton);
0584 
0585       //fill photon pfclusteriso valuemap vectors
0586       int subindex = 0;
0587       for (const auto& ootPhotonFloatValueMapHandle : ootPhotonFloatValueMapHandles) {
0588         ootPhotonFloatValueMapVals[subindex++].push_back((*ootPhotonFloatValueMapHandle)[ootPhotonref]);
0589       }
0590 
0591       //link photon core
0592       const reco::PhotonCoreRef& ootPhotonCore = ootPhoton.photonCore();
0593       linkCore(ootPhotonCore, ootPhotonCores, ootPhotonCoreMap);
0594 
0595       bool slimRelink = slimRelinkOOTPhotonSel_(ootPhoton);
0596       //no supercluster relinking unless slimRelink selection is satisfied
0597       if (!slimRelink)
0598         continue;
0599 
0600       bool relink = relinkOOTPhotonSel_(ootPhoton);
0601 
0602       const reco::SuperClusterRef& ootSuperCluster = ootPhoton.superCluster();
0603       linkSuperCluster(ootSuperCluster, ootSuperClusterMap, ootSuperClusters, relink, ootSuperClusterFullRelinkMap);
0604       //hcal hits
0605       linkHcalHits(*ootPhoton.superCluster(), *hbheHitHandle, hcalRechitMap);
0606     }
0607   }
0608 
0609   //loop over electrons and fill maps
0610   index = -1;
0611   for (const auto& gsfElectron : *gsfElectronHandle) {
0612     index++;
0613 
0614     reco::GsfElectronRef gsfElectronref(gsfElectronHandle, index);
0615     gsfElectrons.push_back(gsfElectron);
0616     auto& newGsfElectron = gsfElectrons.back();
0617     if ((applyGsfElectronCalibOnData_ && event.isRealData()) || (applyGsfElectronCalibOnMC_ && !event.isRealData())) {
0618       calibrateElectron(newGsfElectron,
0619                         gsfElectronref,
0620                         *gsfElectronCalibEnergyHandle,
0621                         *gsfElectronCalibEnergyErrHandle,
0622                         *gsfElectronCalibEcalEnergyHandle,
0623                         *gsfElectronCalibEcalEnergyErrHandle);
0624     }
0625 
0626     bool keep = keepGsfElectronSel_(newGsfElectron);
0627     if (!keep) {
0628       gsfElectrons.pop_back();
0629       continue;
0630     }
0631 
0632     pfCandIsoPairVecEle.push_back((*gsfElectronPfCandMapHandle)[gsfElectronref]);
0633 
0634     //fill electron id valuemap vectors
0635     int subindex = 0;
0636     for (const auto& gsfElectronIdHandle : gsfElectronIdHandles) {
0637       gsfElectronIdVals[subindex++].push_back((*gsfElectronIdHandle)[gsfElectronref]);
0638     }
0639 
0640     subindex = 0;
0641     for (const auto& gsfElectronFloatValueMapHandle : gsfElectronFloatValueMapHandles) {
0642       gsfElectronFloatValueMapVals[subindex++].push_back((*gsfElectronFloatValueMapHandle)[gsfElectronref]);
0643     }
0644 
0645     const reco::GsfElectronCoreRef& gsfElectronCore = gsfElectron.core();
0646     linkCore(gsfElectronCore, gsfElectronCores, gsfElectronCoreMap);
0647 
0648     const reco::GsfTrackRef& gsfTrack = gsfElectron.gsfTrack();
0649 
0650     // Save the main gsfTrack
0651     if (!gsfTrackMap.count(gsfTrack)) {
0652       gsfTracks.push_back(*gsfTrack);
0653       gsfTrackMap[gsfTrack] = gsfTracks.size() - 1;
0654     }
0655 
0656     // Save additional ambiguous gsf tracks in a map:
0657     for (auto const& ambigGsfTrack : gsfElectron.ambiguousGsfTracks()) {
0658       if (!gsfTrackMap.count(ambigGsfTrack)) {
0659         gsfTracks.push_back(*ambigGsfTrack);
0660         gsfTrackMap[ambigGsfTrack] = gsfTracks.size() - 1;
0661       }
0662     }
0663 
0664     bool slimRelink = slimRelinkGsfElectronSel_(newGsfElectron);
0665     //no supercluster relinking unless slimRelink selection is satisfied
0666     if (!slimRelink)
0667       continue;
0668 
0669     bool relink = relinkGsfElectronSel_(newGsfElectron);
0670 
0671     const reco::SuperClusterRef& superCluster = gsfElectron.superCluster();
0672     linkSuperCluster(superCluster, superClusterMap, superClusters, relink, superClusterFullRelinkMap);
0673 
0674     //conversions only for full relinking
0675     if (!relink)
0676       continue;
0677 
0678     const reco::ConversionRefVector& convrefs = gsfElectron.core()->conversions();
0679     linkConversions(convrefs, conversions, conversionMap);
0680 
0681     //explicitly references conversions
0682     const reco::ConversionRefVector& singleconvrefs = gsfElectron.core()->conversionsOneLeg();
0683     linkConversions(singleconvrefs, singleConversions, singleConversionMap);
0684 
0685     //conversions matched by trackrefs
0686     linkConversionsByTrackRef(conversionHandle, gsfElectron, conversions, conversionMap);
0687 
0688     //single leg conversions matched by trackrefs
0689     linkConversionsByTrackRef(singleConversionHandle, gsfElectron, singleConversions, singleConversionMap);
0690 
0691     //hcal hits
0692     linkHcalHits(*gsfElectron.superCluster(), *hbheHitHandle, hcalRechitMap);
0693   }
0694 
0695   //loop over input SuperClusters
0696   index = -1;
0697   for (const auto& superCluster : *scHandle) {
0698     index++;
0699 
0700     const double superclusPt = superCluster.energy() / std::cosh(superCluster.eta());
0701 
0702     if (superclusPt < scPtMin_)
0703       continue;
0704 
0705     if (std::abs(superCluster.eta()) > scAbsetaMax_)
0706       continue;
0707 
0708     bool relinkSupercluster = superclusPt > relinkSuperclusterPtMin_;
0709 
0710     reco::SuperClusterRef superClusterRef(scHandle, index);
0711     linkSuperCluster(superClusterRef, superClusterMap, superClusters, relinkSupercluster, superClusterFullRelinkMap);
0712   }
0713 
0714   //loop over output SuperClusters and fill maps
0715   index = 0;
0716   for (auto& superCluster : superClusters) {
0717     //link seed cluster no matter what
0718     const reco::CaloClusterPtr& seedCluster = superCluster.seed();
0719     linkCaloCluster(seedCluster, ebeeClusters, ebeeClusterMap);
0720 
0721     //only proceed if superCluster is marked for full relinking
0722     bool fullrelink = superClusterFullRelinkMap.count(index++);
0723     if (!fullrelink) {
0724       //zero detid vector which is anyways not useful without stored rechits
0725       superCluster.clearHitsAndFractions();
0726       continue;
0727     }
0728 
0729     // link calo clusters
0730     linkCaloClusters(superCluster,
0731                      ebeeClusters,
0732                      ebeeClusterMap,
0733                      rechitMap,
0734                      barrelHitHandle,
0735                      endcapHitHandle,
0736                      caloTopology,
0737                      esClusters,
0738                      esClusterMap);
0739 
0740     //conversions matched geometrically
0741     linkConversionsByTrackRef(conversionHandle, superCluster, conversions, conversionMap);
0742 
0743     //single leg conversions matched by trackrefs
0744     linkConversionsByTrackRef(singleConversionHandle, superCluster, singleConversions, singleConversionMap);
0745   }
0746 
0747   //loop over output OOTSuperClusters and fill maps
0748   if (!ootPhotonT_.isUninitialized()) {
0749     index = 0;
0750     for (auto& ootSuperCluster : ootSuperClusters) {
0751       //link seed cluster no matter what
0752       const reco::CaloClusterPtr& ootSeedCluster = ootSuperCluster.seed();
0753       linkCaloCluster(ootSeedCluster, ootEbeeClusters, ootEbeeClusterMap);
0754 
0755       //only proceed if ootSuperCluster is marked for full relinking
0756       bool fullrelink = ootSuperClusterFullRelinkMap.count(index++);
0757       if (!fullrelink) {
0758         //zero detid vector which is anyways not useful without stored rechits
0759         ootSuperCluster.clearHitsAndFractions();
0760         continue;
0761       }
0762 
0763       // link calo clusters
0764       linkCaloClusters(ootSuperCluster,
0765                        ootEbeeClusters,
0766                        ootEbeeClusterMap,
0767                        rechitMap,
0768                        barrelHitHandle,
0769                        endcapHitHandle,
0770                        caloTopology,
0771                        ootEsClusters,
0772                        ootEsClusterMap);
0773     }
0774   }
0775   //now finalize and add to the event collections in "reverse" order
0776 
0777   //rechits (fill output collections of rechits to be stored)
0778   for (const EcalRecHit& rechit : *barrelHitHandle) {
0779     if (rechitMap.count(rechit.detid())) {
0780       ebRecHits.push_back(rechit);
0781     }
0782   }
0783 
0784   for (const EcalRecHit& rechit : *endcapHitHandle) {
0785     if (rechitMap.count(rechit.detid())) {
0786       eeRecHits.push_back(rechit);
0787     }
0788   }
0789 
0790   event.emplace(outEBRecHits_, std::move(ebRecHits));
0791   event.emplace(outEERecHits_, std::move(eeRecHits));
0792 
0793   if (doPreshowerEcalHits_) {
0794     for (const EcalRecHit& rechit : *preshowerHitHandle) {
0795       if (rechitMap.count(rechit.detid())) {
0796         esRecHits.push_back(rechit);
0797       }
0798     }
0799     event.emplace(outESRecHits_, std::move(esRecHits));
0800   }
0801 
0802   for (const HBHERecHit& rechit : *hbheHitHandle) {
0803     if (hcalRechitMap.count(rechit.detid())) {
0804       hbheRecHits.push_back(rechit);
0805     }
0806   }
0807   event.emplace(outHBHERecHits_, std::move(hbheRecHits));
0808 
0809   //CaloClusters
0810   //put calocluster output collections in event and get orphan handles to create ptrs
0811   const auto& outEBEEClusterHandle = event.emplace(outEBEEClusters_, std::move(ebeeClusters));
0812   const auto& outESClusterHandle = event.emplace(outESClusters_, std::move(esClusters));
0813   ;
0814 
0815   //Loop over SuperClusters and relink GEDPhoton + GSFElectron CaloClusters
0816   for (reco::SuperCluster& superCluster : superClusters) {
0817     relinkCaloClusters(superCluster, ebeeClusterMap, esClusterMap, outEBEEClusterHandle, outESClusterHandle);
0818   }
0819 
0820   //OOTCaloClusters
0821   //put ootcalocluster output collections in event and get orphan handles to create ptrs
0822   edm::OrphanHandle<reco::CaloClusterCollection> outOOTEBEEClusterHandle;
0823   edm::OrphanHandle<reco::CaloClusterCollection> outOOTESClusterHandle;
0824   //Loop over OOTSuperClusters and relink OOTPhoton CaloClusters
0825   if (!ootPhotonT_.isUninitialized()) {
0826     outOOTEBEEClusterHandle = event.emplace(outOOTEBEEClusters_, std::move(ootEbeeClusters));
0827     outOOTESClusterHandle = event.emplace(outOOTESClusters_, std::move(ootEsClusters));
0828     for (reco::SuperCluster& ootSuperCluster : ootSuperClusters) {
0829       relinkCaloClusters(
0830           ootSuperCluster, ootEbeeClusterMap, ootEsClusterMap, outOOTEBEEClusterHandle, outOOTESClusterHandle);
0831     }
0832   }
0833   //put superclusters and conversions in the event
0834   const auto& outSuperClusterHandle = event.emplace(outSuperClusters_, std::move(superClusters));
0835   const auto& outConversionHandle = event.emplace(outConversions_, std::move(conversions));
0836   const auto& outSingleConversionHandle = event.emplace(outSingleConversions_, std::move(singleConversions));
0837   const auto& outGsfTrackHandle = event.emplace(outGsfTracks_, std::move(gsfTracks));
0838 
0839   //Loop over PhotonCores and relink GEDPhoton SuperClusters (and conversions)
0840   for (reco::PhotonCore& photonCore : photonCores) {
0841     // superclusters
0842     relinkSuperCluster(photonCore, superClusterMap, outSuperClusterHandle);
0843 
0844     //conversions
0845     const reco::ConversionRefVector& convrefs = photonCore.conversions();
0846     relinkConversions(photonCore, convrefs, conversionMap, outConversionHandle);
0847 
0848     //single leg conversions
0849     const reco::ConversionRefVector& singleconvrefs = photonCore.conversionsOneLeg();
0850     relinkConversions(photonCore, singleconvrefs, singleConversionMap, outSingleConversionHandle);
0851   }
0852 
0853   //Relink GSFElectron SuperClusters and main GSF Tracks
0854   for (reco::GsfElectronCore& gsfElectronCore : gsfElectronCores) {
0855     relinkSuperCluster(gsfElectronCore, superClusterMap, outSuperClusterHandle);
0856     relinkGsfTrack(gsfElectronCore, gsfTrackMap, outGsfTrackHandle);
0857   }
0858 
0859   //put ootsuperclusters in the event
0860   edm::OrphanHandle<reco::SuperClusterCollection> outOOTSuperClusterHandle;
0861   if (!ootPhotonT_.isUninitialized())
0862     outOOTSuperClusterHandle = event.emplace(outOOTSuperClusters_, std::move(ootSuperClusters));
0863 
0864   //Relink OOTPhoton SuperClusters
0865   for (reco::PhotonCore& ootPhotonCore : ootPhotonCores) {
0866     relinkSuperCluster(ootPhotonCore, ootSuperClusterMap, outOOTSuperClusterHandle);
0867   }
0868 
0869   //put photoncores and gsfelectroncores into the event
0870   const auto& outPhotonCoreHandle = event.emplace(outPhotonCores_, std::move(photonCores));
0871   edm::OrphanHandle<reco::PhotonCoreCollection> outOOTPhotonCoreHandle;
0872   if (!ootPhotonT_.isUninitialized())
0873     outOOTPhotonCoreHandle = event.emplace(outOOTPhotonCores_, std::move(ootPhotonCores));
0874   const auto& outgsfElectronCoreHandle = event.emplace(outGsfElectronCores_, std::move(gsfElectronCores));
0875 
0876   //loop over photons, oot photons, and electrons and relink the cores
0877   for (reco::Photon& photon : photons) {
0878     relinkPhotonCore(photon, photonCoreMap, outPhotonCoreHandle);
0879   }
0880 
0881   if (!ootPhotonT_.isUninitialized()) {
0882     for (reco::Photon& ootPhoton : ootPhotons) {
0883       relinkPhotonCore(ootPhoton, ootPhotonCoreMap, outOOTPhotonCoreHandle);
0884     }
0885   }
0886 
0887   for (reco::GsfElectron& gsfElectron : gsfElectrons) {
0888     relinkGsfElectronCore(gsfElectron, gsfElectronCoreMap, outgsfElectronCoreHandle);
0889 
0890     // -----
0891     // Also in this loop let's relink ambiguous tracks
0892     std::vector<reco::GsfTrackRef> ambigTracksInThisElectron;
0893     // Here we loop over the ambiguous tracks and save them in a vector
0894     for (auto const& igsf : gsfElectron.ambiguousGsfTracks()) {
0895       ambigTracksInThisElectron.push_back(igsf);
0896     }
0897 
0898     // Now we need to clear them (they are the refs to original collection):
0899     gsfElectron.clearAmbiguousGsfTracks();
0900 
0901     // And here we add them back, now from a new reduced collection:
0902     for (const auto& it : ambigTracksInThisElectron) {
0903       const auto& gsftkmapped = gsfTrackMap.find(it);
0904 
0905       if (gsftkmapped != gsfTrackMap.end()) {
0906         reco::GsfTrackRef gsftkref(outGsfTrackHandle, gsftkmapped->second);
0907         gsfElectron.addAmbiguousGsfTrack(gsftkref);
0908       } else
0909         throw cms::Exception("There must be a problem with linking and mapping of ambiguous gsf tracks...");
0910     }
0911 
0912     if (gsfElectron.ambiguousGsfTracksSize() > 0)
0913       gsfElectron.setAmbiguous(true);  // Set the flag
0914 
0915     ambigTracksInThisElectron.clear();
0916   }
0917 
0918   //(finally) store the output photon and electron collections
0919   const auto& outPhotonHandle = event.emplace(outPhotons_, std::move(photons));
0920   edm::OrphanHandle<reco::PhotonCollection> outOOTPhotonHandle;
0921   if (!ootPhotonT_.isUninitialized())
0922     outOOTPhotonHandle = event.emplace(outOOTPhotons_, std::move(ootPhotons));
0923   const auto& outGsfElectronHandle = event.emplace(outGsfElectrons_, std::move(gsfElectrons));
0924 
0925   //still need to output relinked valuemaps
0926 
0927   //photon pfcand isolation valuemap
0928   edm::ValueMap<std::vector<reco::PFCandidateRef>>::Filler fillerPhotons(photonPfCandMap);
0929   fillerPhotons.insert(outPhotonHandle, pfCandIsoPairVecPho.begin(), pfCandIsoPairVecPho.end());
0930   fillerPhotons.fill();
0931 
0932   //electron pfcand isolation valuemap
0933   edm::ValueMap<std::vector<reco::PFCandidateRef>>::Filler fillerGsfElectrons(gsfElectronPfCandMap);
0934   fillerGsfElectrons.insert(outGsfElectronHandle, pfCandIsoPairVecEle.begin(), pfCandIsoPairVecEle.end());
0935   fillerGsfElectrons.fill();
0936 
0937   event.emplace(outPhotonPfCandMap_, std::move(photonPfCandMap));
0938   event.emplace(outGsfElectronPfCandMap_, std::move(gsfElectronPfCandMap));
0939 
0940   //photon id value maps
0941   index = 0;
0942   for (auto const& vals : photonIdVals) {
0943     emplaceValueMap(outPhotonHandle, vals, event, outPhotonIds_[index++]);
0944   }
0945 
0946   //electron id value maps
0947   index = 0;
0948   for (auto const& vals : gsfElectronIdVals) {
0949     emplaceValueMap(outGsfElectronHandle, vals, event, outGsfElectronIds_[index++]);
0950   }
0951 
0952   // photon iso value maps
0953   index = 0;
0954   for (auto const& vals : photonFloatValueMapVals) {
0955     emplaceValueMap(outPhotonHandle, vals, event, outPhotonFloatValueMaps_[index++]);
0956   }
0957 
0958   if (!ootPhotonT_.isUninitialized()) {
0959     //oot photon iso value maps
0960     index = 0;
0961     for (auto const& vals : ootPhotonFloatValueMapVals) {
0962       emplaceValueMap(outOOTPhotonHandle, vals, event, outOOTPhotonFloatValueMaps_[index++]);
0963     }
0964   }
0965 
0966   // HI photon iso value maps
0967   if (!recoHIPhotonIsolationMapInputToken_.isUninitialized()) {
0968     emplaceValueMap(outPhotonHandle, recoHIPhotonIsolationMapInputVals, event, recoHIPhotonIsolationMapOutputName_);
0969   }
0970 
0971   //electron iso value maps
0972   index = 0;
0973   for (auto const& vals : gsfElectronFloatValueMapVals) {
0974     emplaceValueMap(outGsfElectronHandle, vals, event, outGsfElectronFloatValueMaps_[index++]);
0975   }
0976 }
0977 
0978 template <typename T, typename U>
0979 void ReducedEGProducer::linkCore(const T& core, U& cores, std::map<T, unsigned int>& coreMap) {
0980   if (!coreMap.count(core)) {
0981     cores.push_back(*core);
0982     coreMap[core] = cores.size() - 1;
0983   }
0984 }
0985 
0986 void ReducedEGProducer::linkSuperCluster(const reco::SuperClusterRef& superCluster,
0987                                          std::map<reco::SuperClusterRef, unsigned int>& superClusterMap,
0988                                          reco::SuperClusterCollection& superClusters,
0989                                          const bool relink,
0990                                          std::unordered_set<unsigned int>& superClusterFullRelinkMap) {
0991   const auto& mappedsc = superClusterMap.find(superCluster);
0992   //get index in output collection in order to keep track whether superCluster
0993   //will be subject to full relinking
0994   unsigned int mappedscidx = 0;
0995   if (mappedsc == superClusterMap.end()) {
0996     superClusters.push_back(*superCluster);
0997     mappedscidx = superClusters.size() - 1;
0998     superClusterMap[superCluster] = mappedscidx;
0999   } else {
1000     mappedscidx = mappedsc->second;
1001   }
1002 
1003   //additionally mark supercluster for full relinking
1004   if (relink)
1005     superClusterFullRelinkMap.insert(mappedscidx);
1006 }
1007 
1008 void ReducedEGProducer::linkConversions(const reco::ConversionRefVector& convrefs,
1009                                         reco::ConversionCollection& conversions,
1010                                         std::map<reco::ConversionRef, unsigned int>& conversionMap) {
1011   for (const auto& convref : convrefs) {
1012     linkConversion(convref, conversions, conversionMap);
1013   }
1014 }
1015 
1016 void ReducedEGProducer::linkConversionsByTrackRef(const edm::Handle<reco::ConversionCollection>& conversionHandle,
1017                                                   const reco::GsfElectron& gsfElectron,
1018                                                   reco::ConversionCollection& conversions,
1019                                                   std::map<reco::ConversionRef, unsigned int>& conversionMap) {
1020   int index = 0;
1021   for (const auto& conversion : *conversionHandle) {
1022     reco::ConversionRef convref(conversionHandle, index++);
1023 
1024     bool matched = ConversionTools::matchesConversion(gsfElectron, conversion, true, true);
1025     if (!matched)
1026       continue;
1027 
1028     linkConversion(convref, conversions, conversionMap);
1029   }
1030 }
1031 
1032 void ReducedEGProducer::linkConversionsByTrackRef(const edm::Handle<reco::ConversionCollection>& conversionHandle,
1033                                                   const reco::SuperCluster& superCluster,
1034                                                   reco::ConversionCollection& conversions,
1035                                                   std::map<reco::ConversionRef, unsigned int>& conversionMap) {
1036   int index = 0;
1037   for (const auto& conversion : *conversionHandle) {
1038     reco::ConversionRef convref(conversionHandle, index++);
1039 
1040     bool matched = ConversionTools::matchesConversion(superCluster, conversion, 0.2);
1041     if (!matched)
1042       continue;
1043 
1044     linkConversion(convref, conversions, conversionMap);
1045   }
1046 }
1047 
1048 void ReducedEGProducer::linkConversion(const reco::ConversionRef& convref,
1049                                        reco::ConversionCollection& conversions,
1050                                        std::map<reco::ConversionRef, unsigned int>& conversionMap) {
1051   if (!conversionMap.count(convref)) {
1052     conversions.push_back(*convref);
1053     conversionMap[convref] = conversions.size() - 1;
1054   }
1055 }
1056 
1057 void ReducedEGProducer::linkCaloCluster(const reco::CaloClusterPtr& caloCluster,
1058                                         reco::CaloClusterCollection& caloClusters,
1059                                         std::map<reco::CaloClusterPtr, unsigned int>& caloClusterMap) {
1060   if (!caloClusterMap.count(caloCluster)) {
1061     caloClusters.push_back(*caloCluster);
1062     caloClusterMap[caloCluster] = caloClusters.size() - 1;
1063   }
1064 }
1065 
1066 void ReducedEGProducer::linkCaloClusters(const reco::SuperCluster& superCluster,
1067                                          reco::CaloClusterCollection& ebeeClusters,
1068                                          std::map<reco::CaloClusterPtr, unsigned int>& ebeeClusterMap,
1069                                          std::unordered_set<DetId>& rechitMap,
1070                                          const edm::Handle<EcalRecHitCollection>& barrelHitHandle,
1071                                          const edm::Handle<EcalRecHitCollection>& endcapHitHandle,
1072                                          CaloTopology const& caloTopology,
1073                                          reco::CaloClusterCollection& esClusters,
1074                                          std::map<reco::CaloClusterPtr, unsigned int>& esClusterMap) {
1075   for (const auto& cluster : superCluster.clusters()) {
1076     linkCaloCluster(cluster, ebeeClusters, ebeeClusterMap);
1077 
1078     for (const auto& hitfrac : cluster->hitsAndFractions()) {
1079       rechitMap.insert(hitfrac.first);
1080     }
1081     //make sure to also take all hits in the 5x5 around the max energy xtal
1082     bool barrel = cluster->hitsAndFractions().front().first.subdetId() == EcalBarrel;
1083     const EcalRecHitCollection* rhcol = barrel ? barrelHitHandle.product() : endcapHitHandle.product();
1084     DetId seed = EcalClusterTools::getMaximum(*cluster, rhcol).first;
1085 
1086     std::vector<DetId> dets5x5 =
1087         caloTopology.getSubdetectorTopology(DetId::Ecal, barrel ? EcalBarrel : EcalEndcap)->getWindow(seed, 5, 5);
1088     for (const auto& detid : dets5x5) {
1089       rechitMap.insert(detid);
1090     }
1091   }
1092   for (const auto& cluster : superCluster.preshowerClusters()) {
1093     linkCaloCluster(cluster, esClusters, esClusterMap);
1094 
1095     for (const auto& hitfrac : cluster->hitsAndFractions()) {
1096       rechitMap.insert(hitfrac.first);
1097     }
1098   }
1099 }
1100 
1101 void ReducedEGProducer::linkHcalHits(const reco::SuperCluster& superClus,
1102                                      const HBHERecHitCollection& recHits,
1103                                      std::unordered_set<DetId>& hcalDetIds) {
1104   hcalHitSel_.addDetIds(superClus, recHits, hcalDetIds);
1105 }
1106 
1107 void ReducedEGProducer::relinkCaloClusters(reco::SuperCluster& superCluster,
1108                                            const std::map<reco::CaloClusterPtr, unsigned int>& ebeeClusterMap,
1109                                            const std::map<reco::CaloClusterPtr, unsigned int>& esClusterMap,
1110                                            const edm::OrphanHandle<reco::CaloClusterCollection>& outEBEEClusterHandle,
1111                                            const edm::OrphanHandle<reco::CaloClusterCollection>& outESClusterHandle) {
1112   //remap seed cluster
1113   const auto& seedmapped = ebeeClusterMap.find(superCluster.seed());
1114   if (seedmapped != ebeeClusterMap.end()) {
1115     //make new ptr
1116     reco::CaloClusterPtr clusptr(outEBEEClusterHandle, seedmapped->second);
1117     superCluster.setSeed(clusptr);
1118   }
1119 
1120   //remap all clusters
1121   reco::CaloClusterPtrVector clusters;
1122   for (const auto& cluster : superCluster.clusters()) {
1123     const auto& clustermapped = ebeeClusterMap.find(cluster);
1124     if (clustermapped != ebeeClusterMap.end()) {
1125       //make new ptr
1126       reco::CaloClusterPtr clusptr(outEBEEClusterHandle, clustermapped->second);
1127       clusters.push_back(clusptr);
1128     } else {
1129       //can only relink if all clusters are being relinked, so if one is missing, then skip the relinking completely
1130       clusters.clear();
1131       break;
1132     }
1133   }
1134   if (!clusters.empty()) {
1135     superCluster.setClusters(clusters);
1136   }
1137 
1138   //remap preshower clusters
1139   reco::CaloClusterPtrVector esclusters;
1140   for (const auto& cluster : superCluster.preshowerClusters()) {
1141     const auto& clustermapped = esClusterMap.find(cluster);
1142     if (clustermapped != esClusterMap.end()) {
1143       //make new ptr
1144       reco::CaloClusterPtr clusptr(outESClusterHandle, clustermapped->second);
1145       esclusters.push_back(clusptr);
1146     } else {
1147       //can only relink if all clusters are being relinked, so if one is missing, then skip the relinking completely
1148       esclusters.clear();
1149       break;
1150     }
1151   }
1152   if (!esclusters.empty()) {
1153     superCluster.setPreshowerClusters(esclusters);
1154   }
1155 }
1156 
1157 template <typename T>
1158 void ReducedEGProducer::relinkSuperCluster(
1159     T& core,
1160     const std::map<reco::SuperClusterRef, unsigned int>& superClusterMap,
1161     const edm::OrphanHandle<reco::SuperClusterCollection>& outSuperClusterHandle) {
1162   const auto& scmapped = superClusterMap.find(core.superCluster());
1163   if (scmapped != superClusterMap.end()) {
1164     //make new ref
1165     reco::SuperClusterRef scref(outSuperClusterHandle, scmapped->second);
1166     core.setSuperCluster(scref);
1167   }
1168 }
1169 
1170 void ReducedEGProducer::relinkGsfTrack(reco::GsfElectronCore& gsfElectronCore,
1171                                        const std::map<reco::GsfTrackRef, unsigned int>& gsfTrackMap,
1172                                        const edm::OrphanHandle<reco::GsfTrackCollection>& outGsfTrackHandle) {
1173   const auto& gsftkmapped = gsfTrackMap.find(gsfElectronCore.gsfTrack());
1174   if (gsftkmapped != gsfTrackMap.end()) {
1175     reco::GsfTrackRef gsftkref(outGsfTrackHandle, gsftkmapped->second);
1176     gsfElectronCore.setGsfTrack(gsftkref);
1177   }
1178 }
1179 
1180 void ReducedEGProducer::relinkConversions(reco::PhotonCore& photonCore,
1181                                           const reco::ConversionRefVector& convrefs,
1182                                           const std::map<reco::ConversionRef, unsigned int>& conversionMap,
1183                                           const edm::OrphanHandle<reco::ConversionCollection>& outConversionHandle) {
1184   reco::ConversionRefVector outconvrefs;
1185   for (const auto& convref : convrefs) {
1186     const auto& convmapped = conversionMap.find(convref);
1187     if (convmapped != conversionMap.end()) {
1188       //make new ref
1189       reco::ConversionRef outref(outConversionHandle, convmapped->second);
1190     } else {
1191       //can only relink if all conversions are being relinked, so if one is missing, then skip the relinking completely
1192       outconvrefs.clear();
1193       break;
1194     }
1195   }
1196   if (!outconvrefs.empty()) {
1197     photonCore.setConversions(outconvrefs);
1198   }
1199 }
1200 
1201 void ReducedEGProducer::relinkPhotonCore(reco::Photon& photon,
1202                                          const std::map<reco::PhotonCoreRef, unsigned int>& photonCoreMap,
1203                                          const edm::OrphanHandle<reco::PhotonCoreCollection>& outPhotonCoreHandle) {
1204   const auto& coremapped = photonCoreMap.find(photon.photonCore());
1205   if (coremapped != photonCoreMap.end()) {
1206     //make new ref
1207     reco::PhotonCoreRef coreref(outPhotonCoreHandle, coremapped->second);
1208     photon.setPhotonCore(coreref);
1209   }
1210 }
1211 
1212 void ReducedEGProducer::relinkGsfElectronCore(
1213     reco::GsfElectron& gsfElectron,
1214     const std::map<reco::GsfElectronCoreRef, unsigned int>& gsfElectronCoreMap,
1215     const edm::OrphanHandle<reco::GsfElectronCoreCollection>& outgsfElectronCoreHandle) {
1216   const auto& coremapped = gsfElectronCoreMap.find(gsfElectron.core());
1217   if (coremapped != gsfElectronCoreMap.end()) {
1218     //make new ref
1219     reco::GsfElectronCoreRef coreref(outgsfElectronCoreHandle, coremapped->second);
1220     gsfElectron.setCore(coreref);
1221   }
1222 }
1223 
1224 void ReducedEGProducer::calibratePhoton(reco::Photon& photon,
1225                                         const reco::PhotonRef& oldPhoRef,
1226                                         const edm::ValueMap<float>& energyMap,
1227                                         const edm::ValueMap<float>& energyErrMap) {
1228   float newEnergy = energyMap[oldPhoRef];
1229   float newEnergyErr = energyErrMap[oldPhoRef];
1230   photon.setCorrectedEnergy(reco::Photon::P4type::regression2, newEnergy, newEnergyErr, true);
1231 }
1232 
1233 void ReducedEGProducer::calibrateElectron(reco::GsfElectron& electron,
1234                                           const reco::GsfElectronRef& oldEleRef,
1235                                           const edm::ValueMap<float>& energyMap,
1236                                           const edm::ValueMap<float>& energyErrMap,
1237                                           const edm::ValueMap<float>& ecalEnergyMap,
1238                                           const edm::ValueMap<float>& ecalEnergyErrMap) {
1239   const float newEnergy = energyMap[oldEleRef];
1240   const float newEnergyErr = energyErrMap[oldEleRef];
1241   const float newEcalEnergy = ecalEnergyMap[oldEleRef];
1242   const float newEcalEnergyErr = ecalEnergyErrMap[oldEleRef];
1243 
1244   //make a copy of this as the setCorrectedEcalEnergy call with modifiy the electrons p4
1245   const math::XYZTLorentzVector oldP4 = electron.p4();
1246   const float corr = newEnergy / oldP4.E();
1247 
1248   electron.setCorrectedEcalEnergy(newEcalEnergy);
1249   electron.setCorrectedEcalEnergyError(newEcalEnergyErr);
1250 
1251   math::XYZTLorentzVector newP4{oldP4.x() * corr, oldP4.y() * corr, oldP4.z() * corr, newEnergy};
1252   electron.correctMomentum(newP4, electron.trackMomentumError(), newEnergyErr);
1253 }