Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:07

0001 #include "PTDRElectronID.h"
0002 
0003 void PTDRElectronID::setup(const edm::ParameterSet& conf) {
0004   // Get all the parameters
0005   //baseSetup(conf);
0006 
0007   quality_ = conf.getParameter<std::string>("electronQuality");
0008 
0009   useEoverPIn_ = conf.getParameter<std::vector<int> >("useEoverPIn");
0010   useDeltaEtaIn_ = conf.getParameter<std::vector<int> >("useDeltaEtaIn");
0011   useDeltaPhiIn_ = conf.getParameter<std::vector<int> >("useDeltaPhiIn");
0012   useHoverE_ = conf.getParameter<std::vector<int> >("useHoverE");
0013   useE9overE25_ = conf.getParameter<std::vector<int> >("useE9overE25");
0014   useEoverPOut_ = conf.getParameter<std::vector<int> >("useEoverPOut");
0015   useDeltaPhiOut_ = conf.getParameter<std::vector<int> >("useDeltaPhiOut");
0016   useInvEMinusInvP_ = conf.getParameter<std::vector<int> >("useInvEMinusInvP");
0017   useBremFraction_ = conf.getParameter<std::vector<int> >("useBremFraction");
0018   useSigmaEtaEta_ = conf.getParameter<std::vector<int> >("useSigmaEtaEta");
0019   useSigmaPhiPhi_ = conf.getParameter<std::vector<int> >("useSigmaPhiPhi");
0020   acceptCracks_ = conf.getParameter<std::vector<int> >("acceptCracks");
0021 
0022   if (quality_ == "tight") {
0023     cuts_ = conf.getParameter<edm::ParameterSet>("tightEleIDCuts");
0024     variables_ = 2;
0025   } else if (quality_ == "medium") {
0026     cuts_ = conf.getParameter<edm::ParameterSet>("mediumEleIDCuts");
0027     variables_ = 1;
0028   } else if (quality_ == "loose") {
0029     cuts_ = conf.getParameter<edm::ParameterSet>("looseEleIDCuts");
0030     variables_ = 0;
0031   } else {
0032     throw cms::Exception("Configuration")
0033         << "Invalid electronQuality parameter in PTDElectronID: must be tight, medium or loose.";
0034   }
0035 }
0036 
0037 double PTDRElectronID::result(const reco::GsfElectron* electron, const edm::Event& e, const edm::EventSetup& es) {
0038   //determine which element of the cut arrays in cfi file to read
0039   //depending on the electron classification
0040   int icut = 0;
0041   int elClass = electron->classification();
0042   if (electron->isEB())  //barrel
0043   {
0044     if (elClass == reco::GsfElectron::GOLDEN)
0045       icut = 0;
0046     if (elClass == reco::GsfElectron::BIGBREM)
0047       icut = 1;
0048     //if (elClass == reco::GsfElectron::NARROW)    icut=2;
0049     if (elClass == reco::GsfElectron::SHOWERING)
0050       icut = 3;
0051     if (elClass == reco::GsfElectron::GAP)
0052       icut = 8;
0053   }
0054   if (electron->isEE())  //endcap
0055   {
0056     if (elClass == reco::GsfElectron::GOLDEN)
0057       icut = 4;
0058     if (elClass == reco::GsfElectron::BIGBREM)
0059       icut = 5;
0060     //if (elClass == reco::GsfElectron::NARROW)    icut=6;
0061     if (elClass == reco::GsfElectron::SHOWERING)
0062       icut = 7;
0063     if (elClass == reco::GsfElectron::GAP)
0064       icut = 8;
0065   }
0066   if (elClass == reco::GsfElectron::UNKNOWN) {
0067     edm::LogError("PTDRElectronID") << "Error: unrecognized electron classification ";
0068     return 1.;
0069   }
0070 
0071   if (acceptCracks_[variables_])
0072     if (elClass == reco::GsfElectron::GAP)
0073       return 1.;
0074 
0075   if (useEoverPIn_[variables_]) {
0076     double value = electron->eSuperClusterOverP();
0077     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPInMax");
0078     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPInMin");
0079     if (value < mincut[icut] || value > maxcut[icut])
0080       return 0.;
0081   }
0082 
0083   if (useDeltaEtaIn_[variables_]) {
0084     double value = electron->deltaEtaSuperClusterTrackAtVtx();
0085     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
0086     if (fabs(value) > maxcut[icut])
0087       return 0.;
0088   }
0089 
0090   if (useDeltaPhiIn_[variables_]) {
0091     double value = electron->deltaPhiSuperClusterTrackAtVtx();
0092     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
0093     if (fabs(value) > maxcut[icut])
0094       return 0.;
0095   }
0096 
0097   if (useHoverE_[variables_]) {
0098     double value = electron->hadronicOverEm();
0099     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
0100     if (value > maxcut[icut])
0101       return 0.;
0102   }
0103 
0104   if (useEoverPOut_[variables_]) {
0105     double value = electron->eSeedClusterOverPout();
0106     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
0107     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
0108     if (value < mincut[icut] || value > maxcut[icut])
0109       return 0.;
0110   }
0111 
0112   if (useDeltaPhiOut_[variables_]) {
0113     double value = electron->deltaPhiSeedClusterTrackAtCalo();
0114     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiOut");
0115     if (fabs(value) > maxcut[icut])
0116       return 0.;
0117   }
0118 
0119   if (useInvEMinusInvP_[variables_]) {
0120     double value = (1. / electron->caloEnergy()) - (1. / electron->trackMomentumAtVtx().R());
0121     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("invEMinusInvP");
0122     if (value > maxcut[icut])
0123       return 0.;
0124   }
0125 
0126   if (useBremFraction_[variables_]) {
0127     double value = electron->trackMomentumAtVtx().R() - electron->trackMomentumOut().R();
0128     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("bremFraction");
0129     if (value < mincut[icut])
0130       return 0.;
0131   }
0132 
0133   //EcalClusterLazyTools lazyTools = getClusterShape(e,es);
0134   //std::vector<float> vCov = lazyTools.localCovariances(*(electron->superCluster()->seed())) ;
0135   //std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
0136 
0137   if (useE9overE25_[variables_]) {
0138     double value = electron->r9() * electron->superCluster()->energy() / electron->e5x5();
0139     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("E9overE25");
0140     if (fabs(value) < mincut[icut])
0141       return 0.;
0142   }
0143 
0144   if (useSigmaEtaEta_[variables_]) {
0145     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMax");
0146     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMin");
0147     if (electron->sigmaIetaIeta() < mincut[icut] || electron->sigmaIetaIeta() > maxcut[icut])
0148       return 0.;
0149   }
0150 
0151   if (useSigmaPhiPhi_[variables_]) {
0152     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMin");
0153     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMax");
0154     if (electron->sigmaIphiIphi() < mincut[icut] || electron->sigmaIphiIphi() > maxcut[icut])
0155       return 0.;
0156   }
0157 
0158   return 1.;
0159 }