CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/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 
00086 
00087     indexNConstituents_ = index_type (&bits_, "nConstituents");
00088     indexNEF_ = index_type (&bits_, "NEF");
00089     indexNHF_ = index_type (&bits_, "NHF");
00090     indexCEF_ = index_type (&bits_, "CEF");
00091     indexCHF_ = index_type (&bits_, "CHF");
00092     indexNCH_ = index_type (&bits_, "NCH");
00093 
00094     retInternal_ = getBitTemplate();
00095     
00096   }
00097 
00098 
00099  PFJetIDSelectionFunctor( Version_t version,
00100                           Quality_t quality ) :
00101   version_(version), quality_(quality)
00102  {
00103 
00104     push_back("CHF" );
00105     push_back("NHF" );
00106     push_back("CEF" );
00107     push_back("NEF" );
00108     push_back("NCH" );
00109     push_back("nConstituents");
00110 
00111 
00112     // Set some default cuts for LOOSE, TIGHT
00113     if ( quality_ == LOOSE ) {
00114       set("CHF", 0.0);
00115       set("NHF", 0.99);
00116       set("CEF", 0.99);
00117       set("NEF", 0.99);
00118       set("NCH", 0);
00119       set("nConstituents", 1);
00120     } else if ( quality_ == TIGHT ) {
00121       set("CHF", 0.0);
00122       set("NHF", 0.9);
00123       set("CEF", 0.99);
00124       set("NEF", 0.9);
00125       set("NCH", 0);
00126       set("nConstituents", 1);      
00127     }
00128 
00129 
00130     indexNConstituents_ = index_type (&bits_, "nConstituents");
00131     indexNEF_ = index_type (&bits_, "NEF");
00132     indexNHF_ = index_type (&bits_, "NHF");
00133     indexCEF_ = index_type (&bits_, "CEF");
00134     indexCHF_ = index_type (&bits_, "CHF");
00135     indexNCH_ = index_type (&bits_, "NCH");
00136 
00137     retInternal_ = getBitTemplate();   
00138  }
00139                            
00140 
00141   // 
00142   // Accessor from PAT jets
00143   // 
00144   bool operator()( const pat::Jet & jet, pat::strbitset & ret )  
00145   {
00146     if ( version_ == FIRSTDATA ) {
00147       if ( jet.currentJECLevel() == "Uncorrected" ) 
00148         return firstDataCuts( jet, ret );
00149       else 
00150         return firstDataCuts( jet.correctedJet("Uncorrected"), ret );
00151     }
00152     else {
00153       return false;
00154     }
00155   }
00156   using Selector<pat::Jet>::operator();
00157 
00158   // 
00159   // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID. 
00160   // This can be used with reco quantities. 
00161   // 
00162   bool operator()( reco::PFJet const & jet, 
00163                    pat::strbitset & ret )  
00164   {
00165     if ( version_ == FIRSTDATA ) return firstDataCuts( jet, ret );
00166     else {
00167       return false;
00168     }
00169   }
00170   
00171   // 
00172   // cuts based on craft 08 analysis. 
00173   // 
00174   bool firstDataCuts( reco::Jet const & jet,
00175                       pat::strbitset & ret) 
00176   {    
00177 
00178     ret.set(false);
00179 
00180     // cache some variables
00181     double chf = 0.0;
00182     double nhf = 0.0;
00183     double cef = 0.0;
00184     double nef = 0.0;
00185     int    nch = 0;
00186     int    nconstituents = 0;
00187 
00188     // Have to do this because pat::Jet inherits from reco::Jet but not reco::PFJet
00189     reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>(&jet);
00190     pat::Jet const * patJet = dynamic_cast<pat::Jet const *>(&jet);
00191     reco::BasicJet const * basicJet = dynamic_cast<reco::BasicJet const *>(&jet);
00192 
00193     if ( patJet != 0 ) {
00194 
00195       if ( patJet->isPFJet() ) {
00196         chf = patJet->chargedHadronEnergyFraction();
00197         nhf = ( patJet->neutralHadronEnergy() + patJet->HFHadronEnergy() ) / patJet->energy();
00198         cef = patJet->chargedEmEnergyFraction();
00199         nef = patJet->neutralEmEnergyFraction();
00200         nch = patJet->chargedMultiplicity();
00201         nconstituents = patJet->numberOfDaughters();
00202       } 
00203       // Handle the special case where this is a composed jet for
00204       // subjet analyses
00205       else if ( patJet->isBasicJet() ) {
00206         double e_chf = 0.0;
00207         double e_nhf = 0.0;
00208         double e_cef = 0.0;
00209         double e_nef = 0.0;
00210         nch = 0;
00211         nconstituents = 0;
00212 
00213         for ( reco::Jet::const_iterator ibegin = patJet->begin(),
00214                 iend = patJet->end(), isub = ibegin;
00215               isub != iend; ++isub ) {
00216           reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
00217           e_chf += pfsub->chargedHadronEnergy();
00218           e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
00219           e_cef += pfsub->chargedEmEnergy();
00220           e_nef += pfsub->neutralEmEnergy();
00221           nch += pfsub->chargedMultiplicity();
00222           nconstituents += pfsub->numberOfDaughters();
00223         }
00224         double e = patJet->energy();
00225         if ( e > 0.000001 ) {
00226           chf = e_chf / e;
00227           nhf = e_nhf / e;
00228           cef = e_cef / e;
00229           nef = e_nef / e;
00230         } else {
00231           chf = nhf = cef = nef = 0.0;
00232         }
00233       }
00234     } // end if pat jet
00235     else if ( pfJet != 0 ) {
00236       chf = pfJet->chargedHadronEnergyFraction();
00237       nhf = ( pfJet->neutralHadronEnergy() + pfJet->HFHadronEnergy() ) / pfJet->energy();
00238       cef = pfJet->chargedEmEnergyFraction();
00239       nef = pfJet->neutralEmEnergyFraction();
00240       nch = pfJet->chargedMultiplicity();
00241       nconstituents = pfJet->numberOfDaughters();
00242     } // end if PF jet
00243     // Handle the special case where this is a composed jet for
00244     // subjet analyses
00245     else if ( basicJet != 0 ) {
00246       double e_chf = 0.0;
00247       double e_nhf = 0.0;
00248       double e_cef = 0.0;
00249       double e_nef = 0.0;
00250       nch = 0;
00251       nconstituents = 0;
00252       
00253       for ( reco::Jet::const_iterator ibegin = basicJet->begin(),
00254               iend = patJet->end(), isub = ibegin;
00255             isub != iend; ++isub ) {
00256         reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
00257         e_chf += pfsub->chargedHadronEnergy();
00258         e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
00259         e_cef += pfsub->chargedEmEnergy();
00260         e_nef += pfsub->neutralEmEnergy();
00261         nch += pfsub->chargedMultiplicity();
00262         nconstituents += pfsub->numberOfDaughters();
00263       }
00264       double e = basicJet->energy();
00265       if ( e > 0.000001 ) {
00266         chf = e_chf / e;
00267         nhf = e_nhf / e;
00268         cef = e_cef / e;
00269         nef = e_nef / e;
00270       }
00271     } // end if basic jet
00272 
00273 
00274     // Cuts for all |eta|:
00275     if ( ignoreCut(indexNConstituents_) || nconstituents > cut(indexNConstituents_, int() ) ) passCut( ret, indexNConstituents_);
00276     if ( ignoreCut(indexNEF_)           || ( nef < cut(indexNEF_, double()) ) ) passCut( ret, indexNEF_);
00277     if ( ignoreCut(indexNHF_)           || ( nhf < cut(indexNHF_, double()) ) ) passCut( ret, indexNHF_);    
00278     // Cuts for |eta| < 2.4:
00279     if ( ignoreCut(indexCEF_)           || ( cef < cut(indexCEF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCEF_);
00280     if ( ignoreCut(indexCHF_)           || ( chf > cut(indexCHF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCHF_);
00281     if ( ignoreCut(indexNCH_)           || ( nch > cut(indexNCH_, int())    || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexNCH_);    
00282 
00283     setIgnored( ret );
00284     return (bool)ret;
00285   }
00286   
00287  private: // member variables
00288   
00289   Version_t version_;
00290   Quality_t quality_;
00291 
00292   index_type indexNConstituents_;
00293   index_type indexNEF_;
00294   index_type indexNHF_;
00295   index_type indexCEF_;
00296   index_type indexCHF_;
00297   index_type indexNCH_;
00298   
00299 };
00300 
00301 #endif