CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/PhysicsTools/SelectorUtils/interface/PFJetIDSelectionFunctor.h

Go to the documentation of this file.
00001 #ifndef PhysicsTools_PatUtils_interface_PFJetIDSelectionFunctor_h
00002 #define PhysicsTools_PatUtils_interface_PFJetIDSelectionFunctor_h
00003 
00004 
00022 #include "DataFormats/PatCandidates/interface/Jet.h"
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 
00025 #include "PhysicsTools/SelectorUtils/interface/Selector.h"
00026 
00027 #include <TMath.h>
00028 class PFJetIDSelectionFunctor : public Selector<pat::Jet>  {
00029 
00030  public: // interface
00031 
00032   enum Version_t { FIRSTDATA, N_VERSIONS };
00033   enum Quality_t { LOOSE, TIGHT, N_QUALITY};
00034 
00035   PFJetIDSelectionFunctor() {}
00036   
00037  PFJetIDSelectionFunctor( edm::ParameterSet const & params ) 
00038  {
00039    std::string versionStr = params.getParameter<std::string>("version");
00040    std::string qualityStr = params.getParameter<std::string>("quality");
00041 
00042    if ( versionStr == "FIRSTDATA" ) 
00043      version_ = FIRSTDATA;
00044    else
00045      version_ = FIRSTDATA;  
00046 
00047    if      ( qualityStr == "LOOSE") quality_ = LOOSE;
00048    else if ( qualityStr == "TIGHT") quality_ = TIGHT;
00049    else quality_ = LOOSE;
00050 
00051     push_back("CHF" );
00052     push_back("NHF" );
00053     push_back("CEF" );
00054     push_back("NEF" );
00055     push_back("NCH" );
00056     push_back("nConstituents");
00057 
00058 
00059     // Set some default cuts for LOOSE, TIGHT
00060     if ( quality_ == LOOSE ) {
00061       set("CHF", 0.0);
00062       set("NHF", 0.99);
00063       set("CEF", 0.99);
00064       set("NEF", 0.99);
00065       set("NCH", 0);
00066       set("nConstituents", 1);
00067     } else if ( quality_ == TIGHT ) {
00068       set("CHF", 0.0);
00069       set("NHF", 0.9);
00070       set("CEF", 0.99);
00071       set("NEF", 0.9);
00072       set("NCH", 0);
00073       set("nConstituents", 1);      
00074     }
00075 
00076 
00077     // Now check the configuration to see if the user changed anything
00078     if ( params.exists("CHF") ) set("CHF", params.getParameter<double>("CHF") );
00079     if ( params.exists("NHF") ) set("NHF", params.getParameter<double>("NHF") );
00080     if ( params.exists("CEF") ) set("CEF", params.getParameter<double>("CEF") );
00081     if ( params.exists("NEF") ) set("NEF", params.getParameter<double>("NEF") );
00082     if ( params.exists("NCH") ) set("NCH", params.getParameter<int>   ("NCH") );
00083     if ( params.exists("nConstuents") ) set("nConstituents", params.getParameter<int> ("nConstituents") );
00084 
00085     if ( params.exists("cutsToIgnore") )
00086       setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
00087 
00088 
00089     indexNConstituents_ = index_type (&bits_, "nConstituents");
00090     indexNEF_ = index_type (&bits_, "NEF");
00091     indexNHF_ = index_type (&bits_, "NHF");
00092     indexCEF_ = index_type (&bits_, "CEF");
00093     indexCHF_ = index_type (&bits_, "CHF");
00094     indexNCH_ = index_type (&bits_, "NCH");
00095 
00096     retInternal_ = getBitTemplate();
00097     
00098   }
00099 
00100 
00101  PFJetIDSelectionFunctor( Version_t version,
00102                           Quality_t quality ) :
00103   version_(version), quality_(quality)
00104  {
00105 
00106     push_back("CHF" );
00107     push_back("NHF" );
00108     push_back("CEF" );
00109     push_back("NEF" );
00110     push_back("NCH" );
00111     push_back("nConstituents");
00112 
00113 
00114     // Set some default cuts for LOOSE, TIGHT
00115     if ( quality_ == LOOSE ) {
00116       set("CHF", 0.0);
00117       set("NHF", 0.99);
00118       set("CEF", 0.99);
00119       set("NEF", 0.99);
00120       set("NCH", 0);
00121       set("nConstituents", 1);
00122     } else if ( quality_ == TIGHT ) {
00123       set("CHF", 0.0);
00124       set("NHF", 0.9);
00125       set("CEF", 0.99);
00126       set("NEF", 0.9);
00127       set("NCH", 0);
00128       set("nConstituents", 1);      
00129     }
00130 
00131 
00132     indexNConstituents_ = index_type (&bits_, "nConstituents");
00133     indexNEF_ = index_type (&bits_, "NEF");
00134     indexNHF_ = index_type (&bits_, "NHF");
00135     indexCEF_ = index_type (&bits_, "CEF");
00136     indexCHF_ = index_type (&bits_, "CHF");
00137     indexNCH_ = index_type (&bits_, "NCH");
00138 
00139     retInternal_ = getBitTemplate();   
00140  }
00141                            
00142 
00143   // 
00144   // Accessor from PAT jets
00145   // 
00146   bool operator()( const pat::Jet & jet, pat::strbitset & ret )  
00147   {
00148     if ( version_ == FIRSTDATA ) {
00149       if ( jet.currentJECLevel() == "Uncorrected" ) 
00150         return firstDataCuts( jet, ret );
00151       else 
00152         return firstDataCuts( jet.correctedJet("Uncorrected"), ret );
00153     }
00154     else {
00155       return false;
00156     }
00157   }
00158   using Selector<pat::Jet>::operator();
00159 
00160   // 
00161   // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID. 
00162   // This can be used with reco quantities. 
00163   // 
00164   bool operator()( reco::PFJet const & jet, 
00165                    pat::strbitset & ret )  
00166   {
00167     if ( version_ == FIRSTDATA ) return firstDataCuts( jet, ret );
00168     else {
00169       return false;
00170     }
00171   }
00172   
00173   // 
00174   // cuts based on craft 08 analysis. 
00175   // 
00176   bool firstDataCuts( reco::Jet const & jet,
00177                       pat::strbitset & ret) 
00178   {    
00179 
00180     ret.set(false);
00181 
00182     // cache some variables
00183     double chf = 0.0;
00184     double nhf = 0.0;
00185     double cef = 0.0;
00186     double nef = 0.0;
00187     int    nch = 0;
00188     int    nconstituents = 0;
00189 
00190     // Have to do this because pat::Jet inherits from reco::Jet but not reco::PFJet
00191     reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>(&jet);
00192     pat::Jet const * patJet = dynamic_cast<pat::Jet const *>(&jet);
00193     reco::BasicJet const * basicJet = dynamic_cast<reco::BasicJet const *>(&jet);
00194 
00195     if ( patJet != 0 ) {
00196 
00197       if ( patJet->isPFJet() ) {
00198         chf = patJet->chargedHadronEnergyFraction();
00199         nhf = ( patJet->neutralHadronEnergy() + patJet->HFHadronEnergy() ) / patJet->energy();
00200         cef = patJet->chargedEmEnergyFraction();
00201         nef = patJet->neutralEmEnergyFraction();
00202         nch = patJet->chargedMultiplicity();
00203         nconstituents = patJet->numberOfDaughters();
00204       } 
00205       // Handle the special case where this is a composed jet for
00206       // subjet analyses
00207       else if ( patJet->isBasicJet() ) {
00208         double e_chf = 0.0;
00209         double e_nhf = 0.0;
00210         double e_cef = 0.0;
00211         double e_nef = 0.0;
00212         nch = 0;
00213         nconstituents = 0;
00214 
00215         for ( reco::Jet::const_iterator ibegin = patJet->begin(),
00216                 iend = patJet->end(), isub = ibegin;
00217               isub != iend; ++isub ) {
00218           reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
00219           e_chf += pfsub->chargedHadronEnergy();
00220           e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
00221           e_cef += pfsub->chargedEmEnergy();
00222           e_nef += pfsub->neutralEmEnergy();
00223           nch += pfsub->chargedMultiplicity();
00224           nconstituents += pfsub->numberOfDaughters();
00225         }
00226         double e = patJet->energy();
00227         if ( e > 0.000001 ) {
00228           chf = e_chf / e;
00229           nhf = e_nhf / e;
00230           cef = e_cef / e;
00231           nef = e_nef / e;
00232         } else {
00233           chf = nhf = cef = nef = 0.0;
00234         }
00235       }
00236     } // end if pat jet
00237     else if ( pfJet != 0 ) {
00238       chf = pfJet->chargedHadronEnergyFraction();
00239       nhf = ( pfJet->neutralHadronEnergy() + pfJet->HFHadronEnergy() ) / pfJet->energy();
00240       cef = pfJet->chargedEmEnergyFraction();
00241       nef = pfJet->neutralEmEnergyFraction();
00242       nch = pfJet->chargedMultiplicity();
00243       nconstituents = pfJet->numberOfDaughters();
00244     } // end if PF jet
00245     // Handle the special case where this is a composed jet for
00246     // subjet analyses
00247     else if ( basicJet != 0 ) {
00248       double e_chf = 0.0;
00249       double e_nhf = 0.0;
00250       double e_cef = 0.0;
00251       double e_nef = 0.0;
00252       nch = 0;
00253       nconstituents = 0;
00254       
00255       for ( reco::Jet::const_iterator ibegin = basicJet->begin(),
00256               iend = patJet->end(), isub = ibegin;
00257             isub != iend; ++isub ) {
00258         reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
00259         e_chf += pfsub->chargedHadronEnergy();
00260         e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
00261         e_cef += pfsub->chargedEmEnergy();
00262         e_nef += pfsub->neutralEmEnergy();
00263         nch += pfsub->chargedMultiplicity();
00264         nconstituents += pfsub->numberOfDaughters();
00265       }
00266       double e = basicJet->energy();
00267       if ( e > 0.000001 ) {
00268         chf = e_chf / e;
00269         nhf = e_nhf / e;
00270         cef = e_cef / e;
00271         nef = e_nef / e;
00272       }
00273     } // end if basic jet
00274 
00275 
00276     // Cuts for all |eta|:
00277     if ( ignoreCut(indexNConstituents_) || nconstituents > cut(indexNConstituents_, int() ) ) passCut( ret, indexNConstituents_);
00278     if ( ignoreCut(indexNEF_)           || ( nef < cut(indexNEF_, double()) ) ) passCut( ret, indexNEF_);
00279     if ( ignoreCut(indexNHF_)           || ( nhf < cut(indexNHF_, double()) ) ) passCut( ret, indexNHF_);    
00280     // Cuts for |eta| < 2.4:
00281     if ( ignoreCut(indexCEF_)           || ( cef < cut(indexCEF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCEF_);
00282     if ( ignoreCut(indexCHF_)           || ( chf > cut(indexCHF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCHF_);
00283     if ( ignoreCut(indexNCH_)           || ( nch > cut(indexNCH_, int())    || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexNCH_);    
00284 
00285     setIgnored( ret );
00286     return (bool)ret;
00287   }
00288   
00289  private: // member variables
00290   
00291   Version_t version_;
00292   Quality_t quality_;
00293 
00294   index_type indexNConstituents_;
00295   index_type indexNEF_;
00296   index_type indexNHF_;
00297   index_type indexCEF_;
00298   index_type indexCHF_;
00299   index_type indexNCH_;
00300   
00301 };
00302 
00303 #endif