Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:09

0001 #ifndef PhysicsTools_PatUtils_interface_ElectronVPlusJetsIDSelectionFunctor_h
0002 #define PhysicsTools_PatUtils_interface_ElectronVPlusJetsIDSelectionFunctor_h
0003 
0004 #ifndef __GCCXML__
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 #endif
0007 #include "DataFormats/PatCandidates/interface/Electron.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 
0010 #include "PhysicsTools/SelectorUtils/interface/Selector.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 
0013 class ElectronVPlusJetsIDSelectionFunctor : public Selector<pat::Electron> {
0014 public:  // interface
0015   enum Version_t { SUMMER08, FIRSTDATA, N_VERSIONS };
0016 
0017   ElectronVPlusJetsIDSelectionFunctor() {}
0018 
0019 #ifndef __GCCXML__
0020   ElectronVPlusJetsIDSelectionFunctor(edm::ParameterSet const& parameters, edm::ConsumesCollector& iC)
0021       : ElectronVPlusJetsIDSelectionFunctor(parameters) {}
0022 #endif
0023 
0024   ElectronVPlusJetsIDSelectionFunctor(edm::ParameterSet const& parameters) {
0025     std::string versionStr = parameters.getParameter<std::string>("version");
0026     if (versionStr != "FIRSTDATA") {
0027       std::cout << "The version " << versionStr << " is deprecated. Setting to FIRSTDATA" << std::endl;
0028     }
0029 
0030     if (versionStr == "FIRSTDATA") {
0031       initialize(FIRSTDATA,
0032                  parameters.getParameter<double>("D0"),
0033                  parameters.getParameter<double>("ED0"),
0034                  parameters.getParameter<double>("SD0"),
0035                  parameters.getParameter<double>("RelIso"));
0036       if (parameters.exists("cutsToIgnore"))
0037         setIgnoredCuts(parameters.getParameter<std::vector<std::string> >("cutsToIgnore"));
0038     } else {
0039       throw cms::Exception("InvalidInput") << "Expect version to be one of SUMMER08, FIRSTDATA," << std::endl;
0040     }
0041 
0042     retInternal_ = getBitTemplate();
0043   }
0044 
0045   ElectronVPlusJetsIDSelectionFunctor(
0046       Version_t version, double d0 = 0.2, double ed0 = 999.0, double sd0 = 999.0, double reliso = 0.1) {
0047     initialize(version, d0, ed0, sd0, reliso);
0048   }
0049 
0050   void initialize(Version_t version, double d0, double ed0, double sd0, double reliso) {
0051     version_ = version;
0052 
0053     push_back("D0", d0);
0054     push_back("ED0", ed0);
0055     push_back("SD0", sd0);
0056     push_back("RelIso", reliso);
0057 
0058     // all on by default
0059     set("D0");
0060     set("ED0");
0061     set("SD0");
0062     set("RelIso");
0063 
0064     indexD0_ = index_type(&bits_, "D0");
0065     indexED0_ = index_type(&bits_, "ED0");
0066     indexSD0_ = index_type(&bits_, "SD0");
0067     indexRelIso_ = index_type(&bits_, "RelIso");
0068   }
0069 
0070   // Allow for multiple definitions of the cuts.
0071   bool operator()(const pat::Electron& electron, pat::strbitset& ret) override {
0072     if (version_ == FIRSTDATA)
0073       return firstDataCuts(electron, ret);
0074     else {
0075       return false;
0076     }
0077   }
0078 
0079   using Selector<pat::Electron>::operator();
0080 
0081   // cuts based on craft 08 analysis.
0082   bool firstDataCuts(const pat::Electron& electron, pat::strbitset& ret) {
0083     ret.set(false);
0084     double corr_d0 = electron.dB();
0085     double corr_ed0 = electron.edB();
0086     double corr_sd0 = (corr_ed0 > 0.000000001) ? corr_d0 / corr_ed0 : 999.0;
0087 
0088     double hcalIso = electron.dr03HcalTowerSumEt();
0089     double ecalIso = electron.dr03EcalRecHitSumEt();
0090     double trkIso = electron.dr03TkSumPt();
0091     double et = electron.et();
0092 
0093     double relIso = (ecalIso + hcalIso + trkIso) / et;
0094 
0095     if (fabs(corr_d0) < cut(indexD0_, double()) || ignoreCut(indexD0_))
0096       passCut(ret, indexD0_);
0097     if (fabs(corr_ed0) < cut(indexED0_, double()) || ignoreCut(indexED0_))
0098       passCut(ret, indexED0_);
0099     if (fabs(corr_sd0) < cut(indexSD0_, double()) || ignoreCut(indexSD0_))
0100       passCut(ret, indexSD0_);
0101     if (relIso < cut(indexRelIso_, double()) || ignoreCut(indexRelIso_))
0102       passCut(ret, indexRelIso_);
0103 
0104     setIgnored(ret);
0105     return (bool)ret;
0106   }
0107 
0108 private:  // member variables
0109   Version_t version_;
0110 
0111   index_type indexD0_;
0112   index_type indexED0_;
0113   index_type indexSD0_;
0114   index_type indexRelIso_;
0115 };
0116 
0117 #endif