CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/PhysicsTools/SelectorUtils/interface/PFElectronSelector.h

Go to the documentation of this file.
00001 #ifndef PhysicsTools_PatUtils_interface_PFElectronSelector_h
00002 #define PhysicsTools_PatUtils_interface_PFElectronSelector_h
00003 
00004 #include "DataFormats/PatCandidates/interface/Electron.h"
00005 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 
00008 #include "PhysicsTools/SelectorUtils/interface/Selector.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 
00011 #include <iostream>
00012 
00013 class PFElectronSelector : public Selector<pat::Electron> {
00014   
00015  public: // interface
00016   
00017   bool verbose_;
00018   
00019   enum Version_t { SPRING11, N_VERSIONS };
00020   
00021   PFElectronSelector() {}
00022   
00023   PFElectronSelector( edm::ParameterSet const & parameters ) {
00024     
00025     verbose_ = false;
00026     
00027     std::string versionStr = parameters.getParameter<std::string>("version");
00028     
00029     Version_t version = N_VERSIONS;
00030     
00031     if ( versionStr == "SPRING11" ) {
00032       version = SPRING11;
00033     }
00034     else {
00035       throw cms::Exception("InvalidInput") << "Expect version to be one of SPRING11" << std::endl;
00036     }
00037 
00038 
00039     initialize( version, 
00040                 parameters.getParameter<double>("MVA"),
00041                 parameters.getParameter<double>("D0")  ,
00042                 parameters.getParameter<int>   ("MaxMissingHits"),
00043                 parameters.getParameter<std::string> ("electronIDused"),
00044                 parameters.getParameter<bool>  ("ConversionRejection"),
00045                 parameters.getParameter<double>("PFIso")
00046                 );
00047     if ( parameters.exists("cutsToIgnore") )
00048       setIgnoredCuts( parameters.getParameter<std::vector<std::string> >("cutsToIgnore") );
00049     
00050     retInternal_ = getBitTemplate();
00051 
00052   }
00053   
00054   void initialize( Version_t version,
00055                    double mva = 0.4,
00056                    double d0 = 0.02,
00057                    int nMissingHits = 1,
00058                    std::string eidUsed = "eidTightMC",
00059                    bool convRej = true,
00060                    double pfiso = 0.15 )
00061   {
00062     version_ = version; 
00063     
00064     //    size_t found;
00065     //    found = eidUsed.find("NONE");
00066     //  if ( found != string::npos)
00067     electronIDvalue_ = eidUsed;
00068 
00069     push_back("D0",        d0     );
00070     push_back("MaxMissingHits", nMissingHits  );
00071     push_back("electronID");
00072     push_back("ConversionRejection" );
00073     push_back("PFIso",    pfiso );
00074     push_back("MVA",       mva   );
00075 
00076     set("D0");
00077     set("MaxMissingHits");
00078     set("electronID");
00079     set("ConversionRejection", convRej);
00080     set("PFIso");   
00081     set("MVA");
00082 
00083     indexD0_             = index_type(&bits_, "D0"           );
00084     indexMaxMissingHits_ = index_type(&bits_, "MaxMissingHits" );
00085     indexElectronId_     = index_type(&bits_, "electronID" );
00086     indexConvRej_        = index_type(&bits_, "ConversionRejection" );
00087     indexPFIso_          = index_type(&bits_, "PFIso"       );
00088     indexMVA_            = index_type(&bits_, "MVA"         );
00089   }
00090 
00091   // Allow for multiple definitions of the cuts. 
00092   bool operator()( const pat::Electron & electron, pat::strbitset & ret ) 
00093   { 
00094     if (version_ == SPRING11 ) return spring11Cuts(electron, ret);
00095     else {
00096       return false;
00097     }
00098   }
00099 
00100   using Selector<pat::Electron>::operator();
00101   
00102   // cuts based on top group L+J synchronization exercise
00103   bool spring11Cuts( const pat::Electron & electron, pat::strbitset & ret)
00104   {
00105 
00106     ret.set(false);
00107 
00108     double mva = electron.mva();
00109     double missingHits = electron.gsfTrack()->trackerExpectedHitsInner().numberOfHits() ;
00110     double corr_d0 = electron.dB();
00111 
00112     // in >= 39x conversion rejection variables are accessible from Gsf electron
00113     Double_t dist = electron.convDist(); // default value is -9999 if conversion partner not found
00114     Double_t dcot = electron.convDcot(); // default value is -9999 if conversion partner not found
00115     bool isConv = fabs(dist) < 0.02 && fabs(dcot) < 0.02;
00116 
00117     int bitWiseResults =  (int) electron.electronID( electronIDvalue_ );  
00118     bool electronIDboolean = ((bitWiseResults & 1) == 1 );
00119 
00120     double chIso = electron.userIsolation(pat::PfChargedHadronIso);
00121     double nhIso = electron.userIsolation(pat::PfNeutralHadronIso);
00122     double gIso  = electron.userIsolation(pat::PfGammaIso);
00123     double et    = electron.et() ;
00124 
00125     double pfIso = (chIso + nhIso + gIso) / et;
00126 
00127     if ( missingHits   <= cut(indexMaxMissingHits_,  double()) || ignoreCut(indexMaxMissingHits_)   ) passCut(ret, indexMaxMissingHits_  );
00128     if ( fabs(corr_d0) <  cut(indexD0_,              double()) || ignoreCut(indexD0_)               ) passCut(ret, indexD0_     );
00129     if ( isConv                                                || ignoreCut(indexConvRej_)          ) passCut(ret, indexConvRej_     );
00130     if ( pfIso         <  cut(indexPFIso_,           double()) || ignoreCut(indexPFIso_)            ) passCut(ret, indexPFIso_ );
00131     if ( mva           >  cut(indexMVA_,             double()) || ignoreCut(indexMVA_)              ) passCut(ret, indexMVA_ );
00132     if ( electronIDboolean                                     || ignoreCut(indexElectronId_)       ) passCut(ret, indexElectronId_);
00133     setIgnored(ret);
00134     return (bool)ret;
00135   }
00136   
00137  private: // member variables
00138   
00139   Version_t version_;
00140 
00141   index_type indexID;
00142   index_type indexMaxMissingHits_;
00143   index_type indexD0_;
00144   index_type indexConvRej_;
00145   index_type indexPFIso_;
00146   index_type indexMVA_;
00147   index_type indexElectronId_;
00148 
00149   std::string electronIDvalue_;
00150 };
00151 
00152 #endif