Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:33

0001 //
0002 // Original Authors: Nicholas Wardle, Florian Beaudette
0003 //
0004 #include "RecoParticleFlow/PFProducer/interface/PFEGammaFilters.h"
0005 #include "RecoParticleFlow/PFTracking/interface/PFTrackAlgoTools.h"
0006 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0007 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
0008 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
0009 
0010 using namespace std;
0011 using namespace reco;
0012 
0013 namespace {
0014 
0015   // Constants defining the ECAL barrel limit
0016   constexpr float ecalBarrelMaxEtaWithGap = 1.566;
0017   constexpr float ecalBarrelMaxEtaNoGap = 1.485;
0018 
0019   void readEBEEParams_(const edm::ParameterSet& pset, const std::string& name, std::array<float, 2>& out) {
0020     const auto& vals = pset.getParameter<std::vector<double>>(name);
0021     if (vals.size() != 2)
0022       throw cms::Exception("Configuration") << "Parameter " << name << " does not contain exactly 2 values (EB, EE)\n";
0023     out[0] = vals[0];
0024     out[1] = vals[1];
0025   }
0026 
0027 }  // namespace
0028 
0029 PFEGammaFilters::PFEGammaFilters(const edm::ParameterSet& cfg)
0030     : ph_Et_(cfg.getParameter<double>("photon_MinEt")),
0031       ph_combIso_(cfg.getParameter<double>("photon_combIso")),
0032       ph_loose_hoe_(cfg.getParameter<double>("photon_HoE")),
0033       ph_sietaieta_eb_(cfg.getParameter<double>("photon_SigmaiEtaiEta_barrel")),
0034       ph_sietaieta_ee_(cfg.getParameter<double>("photon_SigmaiEtaiEta_endcap")),
0035       useElePFidDNN_(cfg.getParameter<bool>("useElePFidDnn")),
0036       usePhotonPFidDNN_(cfg.getParameter<bool>("usePhotonPFidDnn")),
0037       useEBModelInGap_(cfg.getParameter<bool>("useEBModelInGap")),
0038       endcapBoundary_(cfg.getParameter<double>("endcapBoundary")),
0039       extEtaBoundary_(cfg.getParameter<double>("extEtaBoundary")),
0040       ele_iso_pt_(cfg.getParameter<double>("electron_iso_pt")),
0041       ele_iso_mva_eb_(cfg.getParameter<double>("electron_iso_mva_barrel")),
0042       ele_iso_mva_ee_(cfg.getParameter<double>("electron_iso_mva_endcap")),
0043       ele_iso_combIso_eb_(cfg.getParameter<double>("electron_iso_combIso_barrel")),
0044       ele_iso_combIso_ee_(cfg.getParameter<double>("electron_iso_combIso_endcap")),
0045       ele_noniso_mva_(cfg.getParameter<double>("electron_noniso_mvaCut")),
0046       ele_missinghits_(cfg.getParameter<unsigned int>("electron_missinghits")),
0047       ele_ecalDrivenHademPreselCut_(cfg.getParameter<double>("electron_ecalDrivenHademPreselCut")),
0048       ele_maxElePtForOnlyMVAPresel_(cfg.getParameter<double>("electron_maxElePtForOnlyMVAPresel")) {
0049   auto const& eleProtectionsForBadHcal = cfg.getParameter<edm::ParameterSet>("electron_protectionsForBadHcal");
0050   auto const& eleProtectionsForJetMET = cfg.getParameter<edm::ParameterSet>("electron_protectionsForJetMET");
0051   auto const& phoProtectionsForBadHcal = cfg.getParameter<edm::ParameterSet>("photon_protectionsForBadHcal");
0052   auto const& phoProtectionsForJetMET = cfg.getParameter<edm::ParameterSet>("photon_protectionsForJetMET");
0053   auto const& eleDNNIdThresholds = cfg.getParameter<edm::ParameterSet>("electronDnnThresholds");
0054   auto const& eleDNNBkgIdThresholds = cfg.getParameter<edm::ParameterSet>("electronDnnBkgThresholds");
0055   auto const& photonDNNIdThresholds = cfg.getParameter<edm::ParameterSet>("photonDnnThresholds");
0056 
0057   pho_sumPtTrackIso_ = phoProtectionsForJetMET.getParameter<double>("sumPtTrackIso");
0058   pho_sumPtTrackIsoSlope_ = phoProtectionsForJetMET.getParameter<double>("sumPtTrackIsoSlope");
0059 
0060   ele_maxNtracks_ = eleProtectionsForJetMET.getParameter<double>("maxNtracks");
0061   ele_maxHcalE_ = eleProtectionsForJetMET.getParameter<double>("maxHcalE");
0062   ele_maxTrackPOverEele_ = eleProtectionsForJetMET.getParameter<double>("maxTrackPOverEele");
0063   ele_maxE_ = eleProtectionsForJetMET.getParameter<double>("maxE");
0064   ele_maxEleHcalEOverEcalE_ = eleProtectionsForJetMET.getParameter<double>("maxEleHcalEOverEcalE");
0065   ele_maxEcalEOverPRes_ = eleProtectionsForJetMET.getParameter<double>("maxEcalEOverPRes");
0066   ele_maxEeleOverPoutRes_ = eleProtectionsForJetMET.getParameter<double>("maxEeleOverPoutRes");
0067   ele_maxHcalEOverP_ = eleProtectionsForJetMET.getParameter<double>("maxHcalEOverP");
0068   ele_maxHcalEOverEcalE_ = eleProtectionsForJetMET.getParameter<double>("maxHcalEOverEcalE");
0069   ele_maxEcalEOverP_1_ = eleProtectionsForJetMET.getParameter<double>("maxEcalEOverP_1");
0070   ele_maxEcalEOverP_2_ = eleProtectionsForJetMET.getParameter<double>("maxEcalEOverP_2");
0071   ele_maxEeleOverPout_ = eleProtectionsForJetMET.getParameter<double>("maxEeleOverPout");
0072   ele_maxDPhiIN_ = eleProtectionsForJetMET.getParameter<double>("maxDPhiIN");
0073 
0074   ele_dnnLowPtThr_ = eleDNNIdThresholds.getParameter<double>("electronDnnLowPtThr");
0075   ele_dnnHighPtBarrelThr_ = eleDNNIdThresholds.getParameter<double>("electronDnnHighPtBarrelThr");
0076   ele_dnnHighPtEndcapThr_ = eleDNNIdThresholds.getParameter<double>("electronDnnHighPtEndcapThr");
0077   ele_dnnExtEta1Thr_ = eleDNNIdThresholds.getParameter<double>("electronDnnExtEta1Thr");
0078   ele_dnnExtEta2Thr_ = eleDNNIdThresholds.getParameter<double>("electronDnnExtEta2Thr");
0079 
0080   ele_dnnBkgLowPtThr_ = eleDNNBkgIdThresholds.getParameter<double>("electronDnnBkgLowPtThr");
0081   ele_dnnBkgHighPtBarrelThr_ = eleDNNBkgIdThresholds.getParameter<double>("electronDnnBkgHighPtBarrelThr");
0082   ele_dnnBkgHighPtEndcapThr_ = eleDNNBkgIdThresholds.getParameter<double>("electronDnnBkgHighPtEndcapThr");
0083   ele_dnnBkgExtEta1Thr_ = eleDNNBkgIdThresholds.getParameter<double>("electronDnnBkgExtEta1Thr");
0084   ele_dnnBkgExtEta2Thr_ = eleDNNBkgIdThresholds.getParameter<double>("electronDnnBkgExtEta2Thr");
0085 
0086   photon_dnnBarrelThr_ = photonDNNIdThresholds.getParameter<double>("photonDnnBarrelThr");
0087   photon_dnnEndcapThr_ = photonDNNIdThresholds.getParameter<double>("photonDnnEndcapThr");
0088 
0089   readEBEEParams_(eleProtectionsForBadHcal, "full5x5_sigmaIetaIeta", badHcal_full5x5_sigmaIetaIeta_);
0090   readEBEEParams_(eleProtectionsForBadHcal, "eInvPInv", badHcal_eInvPInv_);
0091   readEBEEParams_(eleProtectionsForBadHcal, "dEta", badHcal_dEta_);
0092   readEBEEParams_(eleProtectionsForBadHcal, "dPhi", badHcal_dPhi_);
0093   badHcal_eleEnable_ = eleProtectionsForBadHcal.getParameter<bool>("enableProtections");
0094 
0095   badHcal_phoTrkSolidConeIso_offs_ = phoProtectionsForBadHcal.getParameter<double>("solidConeTrkIsoOffset");
0096   badHcal_phoTrkSolidConeIso_slope_ = phoProtectionsForBadHcal.getParameter<double>("solidConeTrkIsoSlope");
0097   badHcal_phoEnable_ = phoProtectionsForBadHcal.getParameter<bool>("enableProtections");
0098 }
0099 
0100 bool PFEGammaFilters::passPhotonSelection(const reco::Photon& photon) const {
0101   // First simple selection, same as the Run1 to be improved in CMSSW_710
0102 
0103   // Photon ET
0104   if (photon.pt() < ph_Et_)
0105     return false;
0106   bool validHoverE = photon.hadTowOverEmValid();
0107   if (debug_)
0108     std::cout << "PFEGammaFilters:: photon pt " << photon.pt() << "   eta, phi " << photon.eta() << ", " << photon.phi()
0109               << "   isoDr03 "
0110               << (photon.trkSumPtHollowConeDR03() + photon.ecalRecHitSumEtConeDR03() + photon.hcalTowerSumEtConeDR03())
0111               << " (cut: " << ph_combIso_ << ")"
0112               << "   H/E " << photon.hadTowOverEm() << " (valid? " << validHoverE << ", cut: " << ph_loose_hoe_ << ")"
0113               << "   s(ieie) " << photon.sigmaIetaIeta()
0114               << " (cut: " << (photon.isEB() ? ph_sietaieta_eb_ : ph_sietaieta_ee_) << ")"
0115               << "   isoTrkDr03Solid " << (photon.trkSumPtSolidConeDR03()) << " (cut: "
0116               << (validHoverE || !badHcal_phoEnable_
0117                       ? -1
0118                       : badHcal_phoTrkSolidConeIso_offs_ + badHcal_phoTrkSolidConeIso_slope_ * photon.pt())
0119               << ")" << std::endl;
0120 
0121   if (usePhotonPFidDNN_) {
0122     // Run3 DNN based PFID
0123     const auto dnn = photon.pfDNN();
0124     const auto photEta = std::abs(photon.eta());
0125     const auto etaThreshold = (useEBModelInGap_) ? ecalBarrelMaxEtaWithGap : ecalBarrelMaxEtaNoGap;
0126     // using the Barrel model for photons in the EB-EE gap
0127     if (photEta <= etaThreshold) {
0128       return dnn > photon_dnnBarrelThr_;
0129     } else if (photEta > etaThreshold) {
0130       return dnn > photon_dnnEndcapThr_;
0131     }
0132   } else {
0133     // Run2 cut based PFID
0134     if (photon.hadTowOverEm() > ph_loose_hoe_)
0135       return false;
0136     //Isolation variables in 0.3 cone combined
0137     if (photon.trkSumPtHollowConeDR03() + photon.ecalRecHitSumEtConeDR03() + photon.hcalTowerSumEtConeDR03() >
0138         ph_combIso_)
0139       return false;
0140 
0141     //patch for bad hcal
0142     if (!validHoverE && badHcal_phoEnable_ &&
0143         photon.trkSumPtSolidConeDR03() >
0144             badHcal_phoTrkSolidConeIso_offs_ + badHcal_phoTrkSolidConeIso_slope_ * photon.pt()) {
0145       return false;
0146     }
0147 
0148     if (photon.isEB()) {
0149       if (photon.sigmaIetaIeta() > ph_sietaieta_eb_)
0150         return false;
0151     } else {
0152       if (photon.sigmaIetaIeta() > ph_sietaieta_ee_)
0153         return false;
0154     }
0155   }
0156 
0157   return true;
0158 }
0159 
0160 bool PFEGammaFilters::passElectronSelection(const reco::GsfElectron& electron,
0161                                             const reco::PFCandidate& pfcand,
0162                                             const int& nVtx) const {
0163   // First simple selection, same as the Run1 to be improved in CMSSW_710
0164 
0165   bool validHoverE = electron.hcalOverEcalValid();
0166   if (debug_)
0167     std::cout << "PFEGammaFilters:: Electron pt " << electron.pt() << " eta, phi " << electron.eta() << ", "
0168               << electron.phi() << " charge " << electron.charge() << " isoDr03 "
0169               << (electron.dr03TkSumPt() + electron.dr03EcalRecHitSumEt() + electron.dr03HcalTowerSumEt())
0170               << " mva_isolated " << electron.mva_Isolated() << " mva_e_pi " << electron.mva_e_pi() << " H/E_valid "
0171               << validHoverE << " s(ieie) " << electron.full5x5_sigmaIetaIeta() << " H/E " << electron.hcalOverEcal()
0172               << " 1/e-1/p " << (1.0 - electron.eSuperClusterOverP()) / electron.ecalEnergy() << " deta "
0173               << std::abs(electron.deltaEtaSeedClusterTrackAtVtx()) << " dphi "
0174               << std::abs(electron.deltaPhiSuperClusterTrackAtVtx()) << endl;
0175 
0176   bool passEleSelection = false;
0177 
0178   // Electron ET
0179   const auto electronPt = electron.pt();
0180   const auto eleEta = std::abs(electron.eta());
0181 
0182   if (useElePFidDNN_) {  // Use DNN for ele pfID >=CMSSW12_1
0183     const auto dnn_sig = electron.dnn_signal_Isolated() + electron.dnn_signal_nonIsolated();
0184     const auto dnn_bkg = electron.dnn_bkg_nonIsolated();
0185     const auto etaThreshold = (useEBModelInGap_) ? ecalBarrelMaxEtaWithGap : ecalBarrelMaxEtaNoGap;
0186     if (eleEta < endcapBoundary_) {
0187       if (electronPt > ele_iso_pt_) {
0188         // using the Barrel model for electron in the EB-EE gap
0189         if (eleEta <= etaThreshold) {  //high pT barrel
0190           passEleSelection = (dnn_sig > ele_dnnHighPtBarrelThr_) && (dnn_bkg < ele_dnnBkgHighPtBarrelThr_);
0191         } else if (eleEta > etaThreshold) {  //high pT endcap (eleEta < 2.5)
0192           passEleSelection = (dnn_sig > ele_dnnHighPtEndcapThr_) && (dnn_bkg < ele_dnnBkgHighPtEndcapThr_);
0193         }
0194       } else {  // pt < ele_iso_pt_ (eleEta < 2.5)
0195         passEleSelection = (dnn_sig > ele_dnnLowPtThr_) && (dnn_bkg < ele_dnnBkgLowPtThr_);
0196       }
0197     } else if ((eleEta >= endcapBoundary_) && (eleEta <= extEtaBoundary_)) {  //First region in extended eta
0198       passEleSelection = (dnn_sig > ele_dnnExtEta1Thr_) && (dnn_bkg < ele_dnnBkgExtEta1Thr_);
0199     } else if (eleEta > extEtaBoundary_) {  //Second region in extended eta
0200       passEleSelection = (dnn_sig > ele_dnnExtEta2Thr_) && (dnn_bkg < ele_dnnBkgExtEta2Thr_);
0201     }
0202     // TODO: For the moment do not evaluate further conditions on isolation and HCAL cleaning..
0203     // To be understood if they are needed
0204   } else {  // Use legacy MVA for ele pfID < CMSSW_12_1
0205     if (electronPt > ele_iso_pt_) {
0206       double isoDr03 = electron.dr03TkSumPt() + electron.dr03EcalRecHitSumEt() + electron.dr03HcalTowerSumEt();
0207       if (eleEta <= ecalBarrelMaxEtaNoGap && isoDr03 < ele_iso_combIso_eb_) {
0208         if (electron.mva_Isolated() > ele_iso_mva_eb_)
0209           passEleSelection = true;
0210       } else if (eleEta > ecalBarrelMaxEtaNoGap && isoDr03 < ele_iso_combIso_ee_) {
0211         if (electron.mva_Isolated() > ele_iso_mva_ee_)
0212           passEleSelection = true;
0213       }
0214     }
0215 
0216     if (electron.mva_e_pi() > ele_noniso_mva_) {
0217       if (validHoverE || !badHcal_eleEnable_) {
0218         passEleSelection = true;
0219       } else {
0220         bool EE = (std::abs(electron.eta()) >
0221                    ecalBarrelMaxEtaNoGap);  // for prefer consistency with above than with E/gamma for now
0222         if ((electron.full5x5_sigmaIetaIeta() < badHcal_full5x5_sigmaIetaIeta_[EE]) &&
0223             (std::abs(1.0 - electron.eSuperClusterOverP()) / electron.ecalEnergy() < badHcal_eInvPInv_[EE]) &&
0224             (std::abs(electron.deltaEtaSeedClusterTrackAtVtx()) <
0225              badHcal_dEta_[EE]) &&  // looser in case of misalignment
0226             (std::abs(electron.deltaPhiSuperClusterTrackAtVtx()) < badHcal_dPhi_[EE])) {
0227           passEleSelection = true;
0228         }
0229       }
0230     }
0231   }
0232   return passEleSelection && passGsfElePreSelWithOnlyConeHadem(electron);
0233 }
0234 
0235 bool PFEGammaFilters::isElectron(const reco::GsfElectron& electron) const {
0236   return electron.gsfTrack()->missingInnerHits() <= ele_missinghits_;
0237 }
0238 
0239 bool PFEGammaFilters::isElectronSafeForJetMET(const reco::GsfElectron& electron,
0240                                               const reco::PFCandidate& pfcand,
0241                                               const reco::Vertex& primaryVertex,
0242                                               bool& lockTracks) const {
0243   bool debugSafeForJetMET = false;
0244   bool isSafeForJetMET = true;
0245 
0246   //   cout << " ele_maxNtracks_ " <<  ele_maxNtracks_ << endl
0247   //        << " ele_maxHcalE_ " << ele_maxHcalE_ << endl
0248   //        << " ele_maxTrackPOverEele_ " << ele_maxTrackPOverEele_ << endl
0249   //        << " ele_maxE_ " << ele_maxE_ << endl
0250   //        << " ele_maxEleHcalEOverEcalE_ "<< ele_maxEleHcalEOverEcalE_ << endl
0251   //        << " ele_maxEcalEOverPRes_ " << ele_maxEcalEOverPRes_ << endl
0252   //        << " ele_maxEeleOverPoutRes_ "  << ele_maxEeleOverPoutRes_ << endl
0253   //        << " ele_maxHcalEOverP_ " << ele_maxHcalEOverP_ << endl
0254   //        << " ele_maxHcalEOverEcalE_ " << ele_maxHcalEOverEcalE_ << endl
0255   //        << " ele_maxEcalEOverP_1_ " << ele_maxEcalEOverP_1_ << endl
0256   //        << " ele_maxEcalEOverP_2_ " << ele_maxEcalEOverP_2_ << endl
0257   //        << " ele_maxEeleOverPout_ "  << ele_maxEeleOverPout_ << endl
0258   //        << " ele_maxDPhiIN_ " << ele_maxDPhiIN_ << endl;
0259 
0260   // loop on the extra-tracks associated to the electron
0261   PFCandidateEGammaExtraRef pfcandextra = pfcand.egammaExtraRef();
0262 
0263   unsigned int iextratrack = 0;
0264   unsigned int itrackHcalLinked = 0;
0265   float SumExtraKfP = 0.;
0266   //float Ene_ecalgsf = 0.;
0267 
0268   // problems here: for now get the electron cluster from the gsf electron
0269   //  const PFCandidate::ElementsInBlocks& eleCluster = pfcandextra->gsfElectronClusterRef();
0270   // PFCandidate::ElementsInBlocks::const_iterator iegfirst = eleCluster.begin();
0271   //  float Ene_hcalgsf = pfcandextra->
0272 
0273   float ETtotal = electron.ecalEnergy();
0274 
0275   //NOTE take this from EGammaExtra
0276   float Ene_ecalgsf = electron.electronCluster()->energy();
0277   float Ene_hcalgsf = pfcandextra->hadEnergy();
0278   float HOverHE = Ene_hcalgsf / (Ene_hcalgsf + Ene_ecalgsf);
0279   float EtotPinMode = electron.eSuperClusterOverP();
0280 
0281   //NOTE take this from EGammaExtra
0282   float EGsfPoutMode = electron.eEleClusterOverPout();
0283   float HOverPin = Ene_hcalgsf / electron.gsfTrack()->pMode();
0284   float dphi_normalsc = electron.deltaPhiSuperClusterTrackAtVtx();
0285 
0286   const PFCandidate::ElementsInBlocks& extraTracks = pfcandextra->extraNonConvTracks();
0287   for (PFCandidate::ElementsInBlocks::const_iterator itrk = extraTracks.begin(); itrk < extraTracks.end(); ++itrk) {
0288     const PFBlock& block = *(itrk->first);
0289     const PFBlock::LinkData& linkData = block.linkData();
0290     const PFBlockElement& pfele = block.elements()[itrk->second];
0291 
0292     if (debugSafeForJetMET)
0293       cout << " My track element number " << itrk->second << endl;
0294     if (pfele.type() == reco::PFBlockElement::TRACK) {
0295       const reco::TrackRef& trackref = pfele.trackRef();
0296 
0297       bool goodTrack = PFTrackAlgoTools::isGoodForEGM(trackref->algo());
0298       // iter0, iter1, iter2, iter3 = Algo < 3
0299       // algo 4,5,6,7
0300 
0301       bool trackIsFromPrimaryVertex = false;
0302       for (Vertex::trackRef_iterator trackIt = primaryVertex.tracks_begin(); trackIt != primaryVertex.tracks_end();
0303            ++trackIt) {
0304         if ((*trackIt).castTo<TrackRef>() == trackref) {
0305           trackIsFromPrimaryVertex = true;
0306           break;
0307         }
0308       }
0309 
0310       // probably we could now remove the algo request??
0311       if (goodTrack && trackref->missingInnerHits() == 0 && trackIsFromPrimaryVertex) {
0312         float p_trk = trackref->p();
0313         SumExtraKfP += p_trk;
0314         iextratrack++;
0315         // Check if these extra tracks are HCAL linked
0316         std::multimap<double, unsigned int> hcalKfElems;
0317         block.associatedElements(
0318             itrk->second, linkData, hcalKfElems, reco::PFBlockElement::HCAL, reco::PFBlock::LINKTEST_ALL);
0319         if (!hcalKfElems.empty()) {
0320           itrackHcalLinked++;
0321         }
0322         if (debugSafeForJetMET)
0323           cout << " The ecalGsf cluster is not isolated: >0 KF extra with algo < 3"
0324                << " Algo " << trackref->algo() << " trackref->missingInnerHits() " << trackref->missingInnerHits()
0325                << " trackIsFromPrimaryVertex " << trackIsFromPrimaryVertex << endl;
0326         if (debugSafeForJetMET)
0327           cout << " My track PT " << trackref->pt() << endl;
0328 
0329       } else {
0330         if (debugSafeForJetMET)
0331           cout << " Tracks from PU "
0332                << " Algo " << trackref->algo() << " trackref->missingInnerHits() " << trackref->missingInnerHits()
0333                << " trackIsFromPrimaryVertex " << trackIsFromPrimaryVertex << endl;
0334         if (debugSafeForJetMET)
0335           cout << " My track PT " << trackref->pt() << endl;
0336       }
0337     }
0338   }
0339   if (iextratrack > 0) {
0340     if (iextratrack > ele_maxNtracks_ || Ene_hcalgsf > ele_maxHcalE_ ||
0341         (SumExtraKfP / Ene_ecalgsf) > ele_maxTrackPOverEele_ ||
0342         (ETtotal > ele_maxE_ && iextratrack > 1 && (Ene_hcalgsf / Ene_ecalgsf) > ele_maxEleHcalEOverEcalE_)) {
0343       if (debugSafeForJetMET)
0344         cout << " *****This electron candidate is discarded: Non isolated  # tracks " << iextratrack << " HOverHE "
0345              << HOverHE << " SumExtraKfP/Ene_ecalgsf " << SumExtraKfP / Ene_ecalgsf << " SumExtraKfP " << SumExtraKfP
0346              << " Ene_ecalgsf " << Ene_ecalgsf << " ETtotal " << ETtotal << " Ene_hcalgsf/Ene_ecalgsf "
0347              << Ene_hcalgsf / Ene_ecalgsf << endl;
0348 
0349       isSafeForJetMET = false;
0350     }
0351     // the electron is retained and the kf tracks are not locked
0352     if ((std::abs(1. - EtotPinMode) < ele_maxEcalEOverPRes_ &&
0353          (std::abs(electron.eta()) < 1.0 || std::abs(electron.eta()) > 2.0)) ||
0354         ((EtotPinMode < 1.1 && EtotPinMode > 0.6) &&
0355          (std::abs(electron.eta()) >= 1.0 && std::abs(electron.eta()) <= 2.0))) {
0356       if (std::abs(1. - EGsfPoutMode) < ele_maxEeleOverPoutRes_ && (itrackHcalLinked == iextratrack)) {
0357         lockTracks = false;
0358         //  lockExtraKf = false;
0359         if (debugSafeForJetMET)
0360           cout << " *****This electron is reactivated  # tracks " << iextratrack << " #tracks hcal linked "
0361                << itrackHcalLinked << " SumExtraKfP/Ene_ecalgsf " << SumExtraKfP / Ene_ecalgsf << " EtotPinMode "
0362                << EtotPinMode << " EGsfPoutMode " << EGsfPoutMode << " eta gsf " << electron.eta() << endl;
0363       }
0364     }
0365   }
0366 
0367   if (HOverPin > ele_maxHcalEOverP_ && HOverHE > ele_maxHcalEOverEcalE_ && EtotPinMode < ele_maxEcalEOverP_1_) {
0368     if (debugSafeForJetMET)
0369       cout << " *****This electron candidate is discarded  HCAL ENERGY "
0370            << " HOverPin " << HOverPin << " HOverHE " << HOverHE << " EtotPinMode" << EtotPinMode << endl;
0371     isSafeForJetMET = false;
0372   }
0373 
0374   // Reject Crazy E/p values... to be understood in the future how to train a
0375   // BDT in order to avoid to select this bad electron candidates.
0376 
0377   if (EtotPinMode < ele_maxEcalEOverP_2_ && EGsfPoutMode < ele_maxEeleOverPout_) {
0378     if (debugSafeForJetMET)
0379       cout << " *****This electron candidate is discarded  Low ETOTPIN "
0380            << " EtotPinMode " << EtotPinMode << " EGsfPoutMode " << EGsfPoutMode << endl;
0381     isSafeForJetMET = false;
0382   }
0383 
0384   // For not-preselected Gsf Tracks ET > 50 GeV, apply dphi preselection
0385   if (ETtotal > ele_maxE_ && std::abs(dphi_normalsc) > ele_maxDPhiIN_) {
0386     if (debugSafeForJetMET)
0387       cout << " *****This electron candidate is discarded  Large ANGLE "
0388            << " ETtotal " << ETtotal << " EGsfPoutMode " << dphi_normalsc << endl;
0389     isSafeForJetMET = false;
0390   }
0391 
0392   return isSafeForJetMET;
0393 }
0394 bool PFEGammaFilters::isPhotonSafeForJetMET(const reco::Photon& photon, const reco::PFCandidate& pfcand) const {
0395   bool isSafeForJetMET = true;
0396   bool debugSafeForJetMET = false;
0397 
0398   //   cout << " pho_sumPtTrackIso_ForPhoton " << pho_sumPtTrackIso_
0399   //        << " pho_sumPtTrackIsoSlope_ForPhoton " << pho_sumPtTrackIsoSlope_ <<  endl;
0400 
0401   float sum_track_pt = 0.;
0402 
0403   PFCandidateEGammaExtraRef pfcandextra = pfcand.egammaExtraRef();
0404   const PFCandidate::ElementsInBlocks& extraTracks = pfcandextra->extraNonConvTracks();
0405   for (PFCandidate::ElementsInBlocks::const_iterator itrk = extraTracks.begin(); itrk < extraTracks.end(); ++itrk) {
0406     const PFBlock& block = *(itrk->first);
0407     const PFBlockElement& pfele = block.elements()[itrk->second];
0408 
0409     if (pfele.type() == reco::PFBlockElement::TRACK) {
0410       const reco::TrackRef& trackref = pfele.trackRef();
0411 
0412       if (debugSafeForJetMET)
0413         cout << "PFEGammaFilters::isPhotonSafeForJetMET photon track:pt " << trackref->pt() << " SingleLegSize "
0414              << pfcandextra->singleLegConvTrackRefMva().size() << endl;
0415 
0416       //const std::vector<reco::TrackRef>&  mySingleLeg =
0417       bool singleLegConv = false;
0418       for (unsigned int iconv = 0; iconv < pfcandextra->singleLegConvTrackRefMva().size(); iconv++) {
0419         if (debugSafeForJetMET)
0420           cout << "PFEGammaFilters::SingleLeg track:pt " << (pfcandextra->singleLegConvTrackRefMva()[iconv].first)->pt()
0421                << endl;
0422 
0423         if (pfcandextra->singleLegConvTrackRefMva()[iconv].first == trackref) {
0424           singleLegConv = true;
0425           if (debugSafeForJetMET)
0426             cout << "PFEGammaFilters::isPhotonSafeForJetMET: SingleLeg conv track " << endl;
0427           break;
0428         }
0429       }
0430       if (singleLegConv)
0431         continue;
0432 
0433       sum_track_pt += trackref->pt();
0434     }
0435   }
0436 
0437   if (debugSafeForJetMET)
0438     cout << " PFEGammaFilters::isPhotonSafeForJetMET: SumPt " << sum_track_pt << endl;
0439 
0440   if (sum_track_pt > (pho_sumPtTrackIso_ + pho_sumPtTrackIsoSlope_ * photon.pt())) {
0441     isSafeForJetMET = false;
0442     if (debugSafeForJetMET)
0443       cout << "************************************!!!! PFEGammaFilters::isPhotonSafeForJetMET: Photon Discaded !!! "
0444            << endl;
0445   }
0446 
0447   return isSafeForJetMET;
0448 }
0449 
0450 //in CMSSW_10_4_0 we changed the electron preselection to be  H/E(cone 0.15) < 0.15
0451 //OR H/E(single tower) < 0.15, with the tower being new.
0452 //However CMS is scared of making any change to the PF content and therefore
0453 //we have to explicitly reject them here
0454 //has to be insync here with GsfElectronAlgo::isPreselected
0455 bool PFEGammaFilters::passGsfElePreSelWithOnlyConeHadem(const reco::GsfElectron& ele) const {
0456   bool passCutBased = ele.passingCutBasedPreselection();
0457   if (ele.hadronicOverEm() > ele_ecalDrivenHademPreselCut_)
0458     passCutBased = false;
0459   bool passMVA = ele.passingMvaPreselection();
0460   if (!ele.ecalDrivenSeed()) {
0461     if (ele.pt() > ele_maxElePtForOnlyMVAPresel_)
0462       return passMVA && passCutBased;
0463     else
0464       return passMVA;
0465   } else
0466     return passCutBased || passMVA;
0467 }
0468 
0469 void PFEGammaFilters::fillPSetDescription(edm::ParameterSetDescription& iDesc) {
0470   // Electron selection cuts
0471   iDesc.add<double>("electron_iso_pt", 10.0);
0472   iDesc.add<double>("electron_iso_mva_barrel", -0.1875);
0473   iDesc.add<double>("electron_iso_mva_endcap", -0.1075);
0474   iDesc.add<double>("electron_iso_combIso_barrel", 10.0);
0475   iDesc.add<double>("electron_iso_combIso_endcap", 10.0);
0476   iDesc.add<double>("electron_noniso_mvaCut", -0.1);
0477   iDesc.add<unsigned int>("electron_missinghits", 1);
0478   iDesc.add<double>("electron_ecalDrivenHademPreselCut", 0.15);
0479   iDesc.add<double>("electron_maxElePtForOnlyMVAPresel", 50.0);
0480   iDesc.add<bool>("useElePFidDnn", false);
0481   iDesc.add<double>("endcapBoundary", 2.5);
0482   iDesc.add<double>("extEtaBoundary", 2.65);
0483   {
0484     edm::ParameterSetDescription psd;
0485     psd.add<double>("electronDnnLowPtThr", 0.5);
0486     psd.add<double>("electronDnnHighPtBarrelThr", 0.5);
0487     psd.add<double>("electronDnnHighPtEndcapThr", 0.5);
0488     psd.add<double>("electronDnnExtEta1Thr", 0.5);
0489     psd.add<double>("electronDnnExtEta2Thr", 0.5);
0490     iDesc.add<edm::ParameterSetDescription>("electronDnnThresholds", psd);
0491   }
0492   {
0493     edm::ParameterSetDescription psd;
0494     psd.add<double>("electronDnnBkgLowPtThr", 1);
0495     psd.add<double>("electronDnnBkgHighPtBarrelThr", 1);
0496     psd.add<double>("electronDnnBkgHighPtEndcapThr", 1);
0497     psd.add<double>("electronDnnBkgExtEta1Thr", 1);
0498     psd.add<double>("electronDnnBkgExtEta2Thr", 1);
0499     iDesc.add<edm::ParameterSetDescription>("electronDnnBkgThresholds", psd);
0500   }
0501   iDesc.add<bool>("usePhotonPFidDnn", false);
0502   {
0503     edm::ParameterSetDescription psd;
0504     psd.add<double>("photonDnnBarrelThr", 0.5);
0505     psd.add<double>("photonDnnEndcapThr", 0.5);
0506     iDesc.add<edm::ParameterSetDescription>("photonDnnThresholds", psd);
0507   }
0508   // control if the EB DNN models should be used up to eta 1.485 or 1.566
0509   iDesc.add<bool>("useEBModelInGap", true);
0510   {
0511     edm::ParameterSetDescription psd;
0512     psd.add<double>("maxNtracks", 3.0)->setComment("Max tracks pointing at Ele cluster");
0513     psd.add<double>("maxHcalE", 10.0);
0514     psd.add<double>("maxTrackPOverEele", 1.0);
0515     psd.add<double>("maxE", 50.0)->setComment("Above this maxE, consider dphi(SC,track) cut");
0516     psd.add<double>("maxEleHcalEOverEcalE", 0.1);
0517     psd.add<double>("maxEcalEOverPRes", 0.2);
0518     psd.add<double>("maxEeleOverPoutRes", 0.5);
0519     psd.add<double>("maxHcalEOverP", 1.0);
0520     psd.add<double>("maxHcalEOverEcalE", 0.1);
0521     psd.add<double>("maxEcalEOverP_1", 0.5)->setComment("E(SC)/P cut - pion rejection");
0522     psd.add<double>("maxEcalEOverP_2", 0.2)->setComment("E(SC)/P cut - weird ele rejection");
0523     psd.add<double>("maxEeleOverPout", 0.2);
0524     psd.add<double>("maxDPhiIN", 0.1)->setComment("Above this dphi(SC,track) and maxE, considered not safe");
0525     iDesc.add<edm::ParameterSetDescription>("electron_protectionsForJetMET", psd);
0526   }
0527   {
0528     edm::ParameterSetDescription psd;
0529     psd.add<bool>("enableProtections", true);
0530     psd.add<std::vector<double>>("full5x5_sigmaIetaIeta",  // EB, EE; 94Xv2 cut-based medium id
0531                                  {0.0106, 0.0387});
0532     psd.add<std::vector<double>>("eInvPInv", {0.184, 0.0721});
0533     psd.add<std::vector<double>>("dEta",  // relax factor 2 to be safer against misalignment
0534                                  {0.0032 * 2, 0.00632 * 2});
0535     psd.add<std::vector<double>>("dPhi", {0.0547, 0.0394});
0536     iDesc.add<edm::ParameterSetDescription>("electron_protectionsForBadHcal", psd);
0537   }
0538 
0539   // Photon selection cuts
0540   iDesc.add<double>("photon_MinEt", 10.0);
0541   iDesc.add<double>("photon_combIso", 10.0);
0542   iDesc.add<double>("photon_HoE", 0.05);
0543   iDesc.add<double>("photon_SigmaiEtaiEta_barrel", 0.0125);
0544   iDesc.add<double>("photon_SigmaiEtaiEta_endcap", 0.034);
0545   {
0546     edm::ParameterSetDescription psd;
0547     psd.add<double>("sumPtTrackIso", 4.0);
0548     psd.add<double>("sumPtTrackIsoSlope", 0.001);
0549     iDesc.add<edm::ParameterSetDescription>("photon_protectionsForJetMET", psd);
0550   }
0551   {
0552     edm::ParameterSetDescription psd;
0553     psd.add<double>("solidConeTrkIsoSlope", 0.3);
0554     psd.add<bool>("enableProtections", true);
0555     psd.add<double>("solidConeTrkIsoOffset", 10.0);
0556     iDesc.add<edm::ParameterSetDescription>("photon_protectionsForBadHcal", psd);
0557   }
0558 }