CMS 3D CMS Logo

Public Member Functions | Protected Attributes

WPlusJetsEventSelector Class Reference

#include <WPlusJetsEventSelector.h>

Inheritance diagram for WPlusJetsEventSelector:
EventSelector

List of all members.

Public Member Functions

std::vector
< reco::ShallowClonePtrCandidate >
const & 
cleanedJets () const
virtual bool operator() (edm::EventBase const &t, pat::strbitset &ret)
void printSelectors (std::ostream &out) const
virtual void scaleJets (double scale)
std::vector
< reco::ShallowClonePtrCandidate >
const & 
selectedElectrons () const
std::vector
< reco::ShallowClonePtrCandidate >
const & 
selectedJets () const
reco::ShallowClonePtrCandidate
const & 
selectedMET () const
std::vector
< reco::ShallowClonePtrCandidate >
const & 
selectedMuons () const
 WPlusJetsEventSelector ()
 WPlusJetsEventSelector (edm::ParameterSet const &params)

Protected Attributes

std::vector
< reco::ShallowClonePtrCandidate
cleanedJets_
index_type conversionIndex_
index_type cosmicIndex_
ElectronVPlusJetsIDSelectionFunctor electronIdLoose_
ElectronVPlusJetsIDSelectionFunctor electronIdTight_
edm::InputTag electronTag_
double eleEtaMax_
double eleEtaMaxLoose_
double eleEtMin_
double eleEtMinLoose_
double eleJetDR_
std::string eleTrig_
bool ePlusJets_
index_type inclusiveIndex_
index_type jet1Index_
index_type jet2Index_
index_type jet3Index_
index_type jet4Index_
index_type jet5Index_
double jetEtaMax_
JetIDSelectionFunctor jetIdLoose_
double jetPtMin_
double jetScale_
edm::InputTag jetTag_
index_type lep1Index_
index_type lep2Index_
index_type lep3Index_
index_type lep4Index_
std::vector
< reco::ShallowClonePtrCandidate
looseElectrons_
std::vector
< reco::ShallowClonePtrCandidate
looseMuons_
reco::ShallowClonePtrCandidate met_
index_type metIndex_
double metMin_
edm::InputTag metTag_
int minJets_
double muEtaMax_
double muEtaMaxLoose_
double muJetDR_
MuonVPlusJetsIDSelectionFunctor muonIdLoose_
MuonVPlusJetsIDSelectionFunctor muonIdTight_
edm::InputTag muonTag_
bool muPlusJets_
double muPtMin_
double muPtMinLoose_
std::string muTrig_
PFJetIDSelectionFunctor pfjetIdLoose_
index_type pvIndex_
PVSelector pvSelector_
std::vector
< reco::ShallowClonePtrCandidate
selectedElectrons2_
std::vector
< reco::ShallowClonePtrCandidate
selectedElectrons_
std::vector
< reco::ShallowClonePtrCandidate
selectedJets_
std::vector
< reco::ShallowClonePtrCandidate
selectedMETs_
std::vector
< reco::ShallowClonePtrCandidate
selectedMuons_
index_type triggerIndex_
edm::InputTag trigTag_
index_type zvetoIndex_

Detailed Description

Definition at line 19 of file WPlusJetsEventSelector.h.


Constructor & Destructor Documentation

WPlusJetsEventSelector::WPlusJetsEventSelector ( ) [inline]

Definition at line 21 of file WPlusJetsEventSelector.h.

{}
WPlusJetsEventSelector::WPlusJetsEventSelector ( edm::ParameterSet const &  params)

Definition at line 9 of file WPlusJetsEventSelector.cc.

References conversionIndex_, cosmicIndex_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), inclusiveIndex_, jet1Index_, jet2Index_, jet3Index_, jet4Index_, jet5Index_, lep1Index_, lep2Index_, lep3Index_, lep4Index_, metIndex_, minJets_, pvIndex_, triggerIndex_, and zvetoIndex_.

                                                                               :
  EventSelector(),
  muonTag_         (params.getParameter<edm::InputTag>("muonSrc") ),
  electronTag_     (params.getParameter<edm::InputTag>("electronSrc") ),
  jetTag_          (params.getParameter<edm::InputTag>("jetSrc") ),
  metTag_          (params.getParameter<edm::InputTag>("metSrc") ),  
  trigTag_         (params.getParameter<edm::InputTag>("trigSrc") ),
  muTrig_          (params.getParameter<std::string>("muTrig")),
  eleTrig_         (params.getParameter<std::string>("eleTrig")),
  pvSelector_      (params.getParameter<edm::ParameterSet>("pvSelector") ),
  muonIdTight_     (params.getParameter<edm::ParameterSet>("muonIdTight") ),
  electronIdTight_ (params.getParameter<edm::ParameterSet>("electronIdTight") ),
  muonIdLoose_     (params.getParameter<edm::ParameterSet>("muonIdLoose") ),
  electronIdLoose_ (params.getParameter<edm::ParameterSet>("electronIdLoose") ),
  jetIdLoose_      (params.getParameter<edm::ParameterSet>("jetIdLoose") ),
  pfjetIdLoose_    (params.getParameter<edm::ParameterSet>("pfjetIdLoose") ),
  minJets_         (params.getParameter<int> ("minJets") ),
  muJetDR_         (params.getParameter<double>("muJetDR")),
  eleJetDR_        (params.getParameter<double>("eleJetDR")),
  muPlusJets_      (params.getParameter<bool>("muPlusJets") ),
  ePlusJets_       (params.getParameter<bool>("ePlusJets") ),
  muPtMin_         (params.getParameter<double>("muPtMin")), 
  muEtaMax_        (params.getParameter<double>("muEtaMax")), 
  eleEtMin_        (params.getParameter<double>("eleEtMin")), 
  eleEtaMax_       (params.getParameter<double>("eleEtaMax")), 
  muPtMinLoose_    (params.getParameter<double>("muPtMinLoose")), 
  muEtaMaxLoose_   (params.getParameter<double>("muEtaMaxLoose")), 
  eleEtMinLoose_   (params.getParameter<double>("eleEtMinLoose")), 
  eleEtaMaxLoose_  (params.getParameter<double>("eleEtaMaxLoose")), 
  jetPtMin_        (params.getParameter<double>("jetPtMin")), 
  jetEtaMax_       (params.getParameter<double>("jetEtaMax")), 
  jetScale_        (params.getParameter<double>("jetScale")),
  metMin_          (params.getParameter<double>("metMin"))
{
  // make the bitset
  push_back( "Inclusive"      );
  push_back( "Trigger"        );
  push_back( "PV"             );
  push_back( ">= 1 Lepton"    );
  push_back( "== 1 Tight Lepton"    );
  push_back( "== 1 Tight Lepton, Mu Veto");
  push_back( "== 1 Lepton"    );
  push_back( "MET Cut"        );
  push_back( "Z Veto"         );
  push_back( "Conversion Veto");
  push_back( "Cosmic Veto"    );
  push_back( ">=1 Jets"       );
  push_back( ">=2 Jets"       );
  push_back( ">=3 Jets"       );
  push_back( ">=4 Jets"       );
  push_back( ">=5 Jets"       );


  // turn (almost) everything on by default
  set( "Inclusive"      );
  set( "Trigger"        );
  set( "PV"             );
  set( ">= 1 Lepton"    );
  set( "== 1 Tight Lepton"    );
  set( "== 1 Tight Lepton, Mu Veto");
  set( "== 1 Lepton"    );
  set( "MET Cut"        );
  set( "Z Veto"         );
  set( "Conversion Veto");
  set( "Cosmic Veto"    );
  set( ">=1 Jets", minJets_ >= 1);
  set( ">=2 Jets", minJets_ >= 2);
  set( ">=3 Jets", minJets_ >= 3);
  set( ">=4 Jets", minJets_ >= 4);
  set( ">=5 Jets", minJets_ >= 5); 


  inclusiveIndex_ = index_type(&bits_, std::string("Inclusive"      ));
  triggerIndex_ = index_type(&bits_, std::string("Trigger"        ));
  pvIndex_ = index_type(&bits_, std::string("PV"             ));
  lep1Index_ = index_type(&bits_, std::string(">= 1 Lepton"    ));
  lep2Index_ = index_type(&bits_, std::string("== 1 Tight Lepton"    ));
  lep3Index_ = index_type(&bits_, std::string("== 1 Tight Lepton, Mu Veto"));
  lep4Index_ = index_type(&bits_, std::string("== 1 Lepton"    ));
  metIndex_ = index_type(&bits_, std::string("MET Cut"        ));
  zvetoIndex_ = index_type(&bits_, std::string("Z Veto"         ));
  conversionIndex_ = index_type(&bits_, std::string("Conversion Veto"));
  cosmicIndex_ = index_type(&bits_, std::string("Cosmic Veto"    ));
  jet1Index_ = index_type(&bits_, std::string(">=1 Jets"));
  jet2Index_ = index_type(&bits_, std::string(">=2 Jets"));
  jet3Index_ = index_type(&bits_, std::string(">=3 Jets"));
  jet4Index_ = index_type(&bits_, std::string(">=4 Jets"));
  jet5Index_ = index_type(&bits_, std::string(">=5 Jets")); 

  if ( params.exists("cutsToIgnore") )
    setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
        

  retInternal_ = getBitTemplate();
}

Member Function Documentation

std::vector<reco::ShallowClonePtrCandidate> const& WPlusJetsEventSelector::cleanedJets ( ) const [inline]

Definition at line 30 of file WPlusJetsEventSelector.h.

References cleanedJets_.

{ return cleanedJets_;      } 
bool WPlusJetsEventSelector::operator() ( edm::EventBase const &  t,
pat::strbitset ret 
) [virtual]

Definition at line 105 of file WPlusJetsEventSelector.cc.

References cleanedJets_, conversionIndex_, cosmicIndex_, reco::deltaR(), electronIdLoose_, electronIdTight_, electronTag_, eleEtaMax_, eleEtaMaxLoose_, eleEtMin_, eleEtMinLoose_, eleJetDR_, eleTrig_, ePlusJets_, reco::LeafCandidate::eta(), Selector< T >::getBitTemplate(), inclusiveIndex_, jet1Index_, jet2Index_, jet3Index_, jet4Index_, jet5Index_, jetEtaMax_, jetIdLoose_, jetPtMin_, jetScale_, jetTag_, lep1Index_, lep2Index_, lep3Index_, lep4Index_, looseElectrons_, looseMuons_, met_, metIndex_, metMin_, metTag_, muEtaMax_, muEtaMaxLoose_, muJetDR_, muonIdLoose_, muonIdTight_, muonTag_, muPlusJets_, muPtMin_, muPtMinLoose_, muTrig_, pat::TriggerEvent::path(), pfjetIdLoose_, reco::LeafCandidate::phi(), reco::LeafCandidate::pt(), pvIndex_, pvSelector_, run_regression::ret, selectedElectrons_, selectedJets_, selectedMETs_, selectedMuons_, pat::strbitset::set(), triggerIndex_, trigTag_, pat::TriggerPath::wasAccept(), pat::TriggerEvent::wasAccept(), pat::TriggerEvent::wasRun(), and zvetoIndex_.

{

  ret.set(false);

  selectedJets_.clear();
  cleanedJets_.clear();
  selectedMuons_.clear();
  selectedElectrons_.clear();
  looseMuons_.clear();
  looseElectrons_.clear();
  selectedMETs_.clear();


  passCut( ret, inclusiveIndex_);


  bool passTrig = false;
  if (!ignoreCut(triggerIndex_) ) {


    edm::Handle<pat::TriggerEvent> triggerEvent;
    event.getByLabel(trigTag_, triggerEvent);

    pat::TriggerEvent const * trig = &*triggerEvent;

    if ( trig->wasRun() && trig->wasAccept() ) {

      pat::TriggerPath const * muPath = trig->path(muTrig_);

      pat::TriggerPath const * elePath = trig->path(eleTrig_);

      if ( muPlusJets_ && muPath != 0 && muPath->wasAccept() ) {
        passTrig = true;    
      }
      
      if ( ePlusJets_ && elePath != 0 && elePath->wasAccept() ) {
        passTrig = true;
      }
    }
  }



  
  if ( ignoreCut(triggerIndex_) || 
       passTrig ) {
    passCut(ret, triggerIndex_);


    bool passPV = false;

    passPV = pvSelector_( event );
    if ( ignoreCut(pvIndex_) || passPV ) {
      passCut(ret, pvIndex_);
  
      edm::Handle< vector< pat::Electron > > electronHandle;
      event.getByLabel (electronTag_, electronHandle);
  
      edm::Handle< vector< pat::Muon > > muonHandle;
      event.getByLabel (muonTag_, muonHandle);

      edm::Handle< vector< pat::Jet > > jetHandle;

      edm::Handle< edm::OwnVector<reco::Candidate> > jetClonesHandle ;

      edm::Handle< vector< pat::MET > > metHandle;
      event.getByLabel (metTag_, metHandle);

      int nElectrons = 0;
      for ( std::vector<pat::Electron>::const_iterator electronBegin = electronHandle->begin(),
              electronEnd = electronHandle->end(), ielectron = electronBegin;
            ielectron != electronEnd; ++ielectron ) {
        ++nElectrons;
        // Tight cuts
        if ( ielectron->et() > eleEtMin_ && fabs(ielectron->eta()) < eleEtaMax_ && 
             electronIdTight_(*ielectron) &&
             ielectron->electronID( "eidRobustTight" ) > 0  ) {
          selectedElectrons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Electron>( electronHandle, ielectron - electronBegin ) ) );
        } else {
          // Loose cuts
          if ( ielectron->et() > eleEtMinLoose_ && fabs(ielectron->eta()) < eleEtaMaxLoose_ && 
               electronIdLoose_(*ielectron) ) {
            looseElectrons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Electron>( electronHandle, ielectron - electronBegin ) ) );
          }
        }
      }


      for ( std::vector<pat::Muon>::const_iterator muonBegin = muonHandle->begin(),
              muonEnd = muonHandle->end(), imuon = muonBegin;
            imuon != muonEnd; ++imuon ) {
        if ( !imuon->isGlobalMuon() ) continue;
        
        // Tight cuts
        bool passTight = muonIdTight_(*imuon,event) && imuon->isTrackerMuon() ;
        if (  imuon->pt() > muPtMin_ && fabs(imuon->eta()) < muEtaMax_ && 
             passTight ) {

          selectedMuons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Muon>( muonHandle, imuon - muonBegin ) ) );
        } else {
          // Loose cuts
          if ( imuon->pt() > muPtMinLoose_ && fabs(imuon->eta()) < muEtaMaxLoose_ && 
               muonIdLoose_(*imuon,event) ) {
            looseMuons_.push_back( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Muon>( muonHandle, imuon - muonBegin ) ) );
          }
        }
      }


      met_ = reco::ShallowClonePtrCandidate( edm::Ptr<pat::MET>( metHandle, 0),
                                             metHandle->at(0).charge(),
                                             metHandle->at(0).p4() );



      event.getByLabel (jetTag_, jetHandle);
      pat::strbitset ret1 = jetIdLoose_.getBitTemplate();
      pat::strbitset ret2 = pfjetIdLoose_.getBitTemplate();
      for ( std::vector<pat::Jet>::const_iterator jetBegin = jetHandle->begin(),
              jetEnd = jetHandle->end(), ijet = jetBegin;
            ijet != jetEnd; ++ijet ) {
        reco::ShallowClonePtrCandidate scaledJet ( reco::ShallowClonePtrCandidate( edm::Ptr<pat::Jet>( jetHandle, ijet - jetBegin ),
                                                                                   ijet->charge(),
                                                                                   ijet->p4() * jetScale_ ) );    
        bool passJetID = false;
        if ( ijet->isCaloJet() || ijet->isJPTJet() ) passJetID = jetIdLoose_(*ijet, ret1);
        else passJetID = pfjetIdLoose_(*ijet, ret2);
        if ( scaledJet.pt() > jetPtMin_ && fabs(scaledJet.eta()) < jetEtaMax_ && passJetID ) {
          selectedJets_.push_back( scaledJet );

          if ( muPlusJets_ ) {

            //Remove some jets
            bool indeltaR = false;
            for( std::vector<reco::ShallowClonePtrCandidate>::const_iterator muonBegin = selectedMuons_.begin(),
                   muonEnd = selectedMuons_.end(), imuon = muonBegin;
                 imuon != muonEnd; ++imuon ) {
              if( reco::deltaR( imuon->eta(), imuon->phi(), scaledJet.eta(), scaledJet.phi() ) < muJetDR_ )
                {  indeltaR = true; }
            }
            if( !indeltaR ) {
              cleanedJets_.push_back( scaledJet );
            }// end if jet is not within dR of a muon
          }// end if mu+jets
          else {
            //Remove some jets
            bool indeltaR = false;
            for( std::vector<reco::ShallowClonePtrCandidate>::const_iterator electronBegin = selectedElectrons_.begin(),
                   electronEnd = selectedElectrons_.end(), ielectron = electronBegin;
                 ielectron != electronEnd; ++ielectron ) {
              if( reco::deltaR( ielectron->eta(), ielectron->phi(), scaledJet.eta(), scaledJet.phi() ) < eleJetDR_ )
                {  indeltaR = true; }
            }
            if( !indeltaR ) {
              cleanedJets_.push_back( scaledJet );
            }// end if jet is not within dR of an electron
          }// end if e+jets
        }// end if pass id and kin cuts
      }// end loop over jets



      int nleptons = 0;
      if ( muPlusJets_ )
        nleptons += selectedMuons_.size();
      
      if ( ePlusJets_ ) 
        nleptons += selectedElectrons_.size();

      if ( ignoreCut(lep1Index_) || 
           ( nleptons > 0 ) ){
        passCut( ret, lep1Index_);

        if ( ignoreCut(lep2Index_) || 
             ( nleptons == 1 ) ){
          passCut( ret, lep2Index_);

          bool oneMuon = 
            ( selectedMuons_.size() == 1 && 
              looseMuons_.size() + selectedElectrons_.size() + looseElectrons_.size() == 0 
              );
          bool oneElectron = 
            ( selectedElectrons_.size() == 1 &&
              selectedMuons_.size() == 0 
              );

          bool oneMuonMuVeto = 
            ( selectedMuons_.size() == 1 && 
              looseMuons_.size() == 0 
              );


          if ( ignoreCut(lep3Index_) || 
               ePlusJets_ ||
               (muPlusJets_ && oneMuonMuVeto)
               ) {
            passCut(ret, lep3Index_);

            if ( ignoreCut(lep4Index_) || 
                 ( (muPlusJets_ && oneMuon) ^ (ePlusJets_ && oneElectron )  )
                 ) {
              passCut(ret, lep4Index_);   

              bool metCut = met_.pt() > metMin_;
              if ( ignoreCut(metIndex_) ||
                   metCut ) {
                passCut( ret, metIndex_ );
          

                bool zVeto = true;
                if ( selectedMuons_.size() == 2 ) {
                }
                if ( selectedElectrons_.size() == 2 ) {
                }
                if ( ignoreCut(zvetoIndex_) ||
                     zVeto ){
                  passCut(ret, zvetoIndex_);
            
  
                  bool conversionVeto = true;
                  if ( ignoreCut(conversionIndex_) ||
                       conversionVeto ) {
                    passCut(ret,conversionIndex_);
                


                    bool cosmicVeto = true;
                    if ( ignoreCut(cosmicIndex_) ||
                         cosmicVeto ) {
                      passCut(ret,cosmicIndex_);

                      if ( ignoreCut(jet1Index_) ||
                           static_cast<int>(cleanedJets_.size()) >=  1 ){
                        passCut(ret,jet1Index_);  
                      } // end if >=1 tight jets

                      if ( ignoreCut(jet2Index_) ||
                           static_cast<int>(cleanedJets_.size()) >=  2 ){
                        passCut(ret,jet2Index_);  
                      } // end if >=2 tight jets

                      if ( ignoreCut(jet3Index_) ||
                           static_cast<int>(cleanedJets_.size()) >=  3 ){
                        passCut(ret,jet3Index_);  
                      } // end if >=3 tight jets

                      if ( ignoreCut(jet4Index_) ||
                           static_cast<int>(cleanedJets_.size()) >=  4 ){
                        passCut(ret,jet4Index_);  
                      } // end if >=4 tight jets

                      if ( ignoreCut(jet5Index_) ||
                           static_cast<int>(cleanedJets_.size()) >=  5 ){
                        passCut(ret,jet5Index_);  
                      } // end if >=5 tight jets


                  
                    } // end if cosmic veto
                
                  } // end if conversion veto

                } // end if z veto

              } // end if met cut
        
            } // end if == 1 lepton

          } // end if == 1 tight lepton with a muon veto separately

        } // end if == 1 tight lepton

      } // end if >= 1 lepton

    } // end if PV
    
  } // end if trigger


  setIgnored(ret);

  return (bool)ret;
}
void WPlusJetsEventSelector::printSelectors ( std::ostream &  out) const [inline]

Definition at line 35 of file WPlusJetsEventSelector.h.

References electronIdLoose_, electronIdTight_, jetIdLoose_, muonIdLoose_, muonIdTight_, pfjetIdLoose_, Selector< T >::print(), and pvSelector_.

                                              {
    out << "PV Selector: " << std::endl;
    pvSelector_.print(out);
    out << "Muon ID Tight Selector: " << std::endl;
    muonIdTight_.print(out);
    out << "Electron ID Tight Selector: " << std::endl;
    electronIdTight_.print(out);
    out << "Muon ID Loose Selector: " << std::endl;
    muonIdLoose_.print(out);
    out << "Electron ID Loose Selector: " << std::endl;
    electronIdLoose_.print(out);
    out << "Calo Jet Selector: " << std::endl;
    jetIdLoose_.print(out);
    out << "PF Jet Selector: " << std::endl;
    pfjetIdLoose_.print(out);
  }
virtual void WPlusJetsEventSelector::scaleJets ( double  scale) [inline, virtual]

Definition at line 24 of file WPlusJetsEventSelector.h.

References jetScale_, and pileupReCalc_HLTpaths::scale.

std::vector<reco::ShallowClonePtrCandidate> const& WPlusJetsEventSelector::selectedElectrons ( ) const [inline]

Definition at line 31 of file WPlusJetsEventSelector.h.

References selectedElectrons_.

{ return selectedElectrons_;}
std::vector<reco::ShallowClonePtrCandidate> const& WPlusJetsEventSelector::selectedJets ( ) const [inline]

Definition at line 29 of file WPlusJetsEventSelector.h.

References selectedJets_.

{ return selectedJets_;     } 
reco::ShallowClonePtrCandidate const& WPlusJetsEventSelector::selectedMET ( ) const [inline]

Definition at line 33 of file WPlusJetsEventSelector.h.

References met_.

{ return met_; }
std::vector<reco::ShallowClonePtrCandidate> const& WPlusJetsEventSelector::selectedMuons ( ) const [inline]

Definition at line 32 of file WPlusJetsEventSelector.h.

References selectedMuons_.

{ return selectedMuons_;    }

Member Data Documentation

Definition at line 69 of file WPlusJetsEventSelector.h.

Referenced by cleanedJets(), and operator()().

Definition at line 115 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::cosmicIndex_ [protected]

Definition at line 116 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

Definition at line 77 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

Definition at line 75 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

Definition at line 55 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 92 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 97 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 91 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 96 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 84 of file WPlusJetsEventSelector.h.

Referenced by operator()().

std::string WPlusJetsEventSelector::eleTrig_ [protected]

Definition at line 61 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 87 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 106 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::jet1Index_ [protected]

Definition at line 117 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::jet2Index_ [protected]

Definition at line 118 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::jet3Index_ [protected]

Definition at line 119 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::jet4Index_ [protected]

Definition at line 120 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::jet5Index_ [protected]

Definition at line 121 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

Definition at line 100 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 78 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

Definition at line 99 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 102 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and scaleJets().

Definition at line 56 of file WPlusJetsEventSelector.h.

Referenced by operator()().

index_type WPlusJetsEventSelector::lep1Index_ [protected]

Definition at line 109 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::lep2Index_ [protected]

Definition at line 110 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::lep3Index_ [protected]

Definition at line 111 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

index_type WPlusJetsEventSelector::lep4Index_ [protected]

Definition at line 112 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

Definition at line 67 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 66 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 71 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and selectedMET().

index_type WPlusJetsEventSelector::metIndex_ [protected]

Definition at line 113 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

double WPlusJetsEventSelector::metMin_ [protected]

Definition at line 104 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 57 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 81 of file WPlusJetsEventSelector.h.

Referenced by WPlusJetsEventSelector().

Definition at line 90 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 95 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 83 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 76 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

Definition at line 74 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

Definition at line 54 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 86 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 89 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 94 of file WPlusJetsEventSelector.h.

Referenced by operator()().

std::string WPlusJetsEventSelector::muTrig_ [protected]

Definition at line 60 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 79 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

index_type WPlusJetsEventSelector::pvIndex_ [protected]

Definition at line 108 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

Definition at line 73 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and printSelectors().

Definition at line 70 of file WPlusJetsEventSelector.h.

Definition at line 65 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and selectedElectrons().

Definition at line 63 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and selectedJets().

Definition at line 68 of file WPlusJetsEventSelector.h.

Referenced by operator()().

Definition at line 64 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and selectedMuons().

index_type WPlusJetsEventSelector::triggerIndex_ [protected]

Definition at line 107 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().

Definition at line 58 of file WPlusJetsEventSelector.h.

Referenced by operator()().

index_type WPlusJetsEventSelector::zvetoIndex_ [protected]

Definition at line 114 of file WPlusJetsEventSelector.h.

Referenced by operator()(), and WPlusJetsEventSelector().