Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "ClassBasedElectronID.h"
0002 
0003 // ===========================================================================================================
0004 void ClassBasedElectronID::setup(const edm::ParameterSet& conf)
0005 // ===========================================================================================================
0006 {
0007   // Get all the parameters
0008   //baseSetup(conf);
0009 
0010   quality_ = conf.getParameter<std::string>("electronQuality");
0011 
0012   if (quality_ == "Eff95Cuts") {
0013     cuts_ = conf.getParameter<edm::ParameterSet>("Eff95Cuts");
0014   }
0015 
0016   else if (quality_ == "Eff90Cuts") {
0017     cuts_ = conf.getParameter<edm::ParameterSet>("Eff90Cuts");
0018   }
0019 
0020   else {
0021     edm::LogError("ClassBasedElectronID") << "Invalid electronQuality parameter: must be tight, medium or loose.";
0022     exit(1);
0023   }
0024 
0025 }  // end of setup
0026 
0027 double ClassBasedElectronID::result(const reco::GsfElectron* electron, const edm::Event& e, const edm::EventSetup& es) {
0028   //determine which element of the cut arrays in cfi file to read
0029   //depending on the electron classification
0030   int icut = 0;
0031   int elClass = electron->classification();
0032   if (electron->isEB())  //barrel
0033   {
0034     if (elClass == reco::GsfElectron::GOLDEN)
0035       icut = 0;
0036     if (elClass == reco::GsfElectron::BIGBREM)
0037       icut = 1;
0038     if (elClass == reco::GsfElectron::SHOWERING)
0039       icut = 2;
0040     if (elClass == reco::GsfElectron::GAP)
0041       icut = 6;
0042   }
0043   if (electron->isEE())  //endcap
0044   {
0045     if (elClass == reco::GsfElectron::GOLDEN)
0046       icut = 3;
0047     if (elClass == reco::GsfElectron::BIGBREM)
0048       icut = 4;
0049     if (elClass == reco::GsfElectron::SHOWERING)
0050       icut = 5;
0051     if (elClass == reco::GsfElectron::GAP)
0052       icut = 7;
0053   }
0054   if (elClass == reco::GsfElectron::UNKNOWN) {
0055     edm::LogError("ClassBasedElectronID") << "Error: unrecognized electron classification ";
0056     return 1.;
0057   }
0058 
0059   bool useDeltaEtaIn = true;
0060   bool useSigmaIetaIeta = true;
0061   bool useHoverE = true;
0062   bool useEoverPOut = true;
0063   bool useDeltaPhiInCharge = true;
0064 
0065   // DeltaEtaIn
0066   if (useDeltaEtaIn) {
0067     double value = electron->deltaEtaSuperClusterTrackAtVtx();
0068     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
0069     if (fabs(value) > maxcut[icut])
0070       return 0.;
0071   }
0072 
0073   // SigmaIetaIeta
0074   if (useSigmaIetaIeta) {
0075     double value = electron->sigmaIetaIeta();
0076     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaIetaIetaMax");
0077     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaIetaIetaMin");
0078     if (value < mincut[icut] || value > maxcut[icut])
0079       return 0.;
0080   }
0081 
0082   // H/E
0083   if (useHoverE) {  //_[variables_]) {
0084     double value = electron->hadronicOverEm();
0085     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
0086     if (value > maxcut[icut])
0087       return 0.;
0088   }  // if use
0089 
0090   // Eseed/Pout
0091   if (useEoverPOut) {
0092     double value = electron->eSeedClusterOverPout();
0093     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
0094     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
0095     if (value < mincut[icut] || value > maxcut[icut])
0096       return 0.;
0097   }
0098 
0099   // DeltaPhiIn*Charge
0100   if (useDeltaPhiInCharge) {
0101     double value1 = electron->deltaPhiSuperClusterTrackAtVtx();
0102     double value2 = electron->charge();
0103     double value = value1 * value2;
0104     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiInChargeMax");
0105     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("deltaPhiInChargeMin");
0106     if (value < mincut[icut] || value > maxcut[icut])
0107       return 0.;
0108   }
0109 
0110   return 1.;
0111 }