CMS 3D CMS Logo

Public Member Functions | Private Attributes

HiggsToZZ4LeptonsSkim Class Reference

#include <HiggsToZZ4LeptonsSkim.h>

Inheritance diagram for HiggsToZZ4LeptonsSkim:
edm::EDFilter edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 Get event properties to send to builder to fill seed collection.
 HiggsToZZ4LeptonsSkim (const edm::ParameterSet &)
 ~HiggsToZZ4LeptonsSkim ()

Private Attributes

bool debug
int nEvents
int nLeptonMin
int nSelectedEvents
int nStiffLeptonMin
edm::InputTag recTrackLabel
float softMinPt
float stiffMinPt
edm::InputTag theGLBMuonLabel
edm::InputTag theGsfELabel

Detailed Description

Definition at line 28 of file HiggsToZZ4LeptonsSkim.h.


Constructor & Destructor Documentation

HiggsToZZ4LeptonsSkim::HiggsToZZ4LeptonsSkim ( const edm::ParameterSet pset) [explicit]

Definition at line 34 of file HiggsToZZ4LeptonsSkim.cc.

References debug, edm::ParameterSet::getParameter(), and nEvents.

                                                                        {

  // Local Debug flag
  debug              = pset.getParameter<bool>("DebugHiggsToZZ4LeptonsSkim");

  // Reconstructed objects
  recTrackLabel      = pset.getParameter<edm::InputTag>("RecoTrackLabel");
  theGLBMuonLabel    = pset.getParameter<edm::InputTag>("GlobalMuonCollectionLabel");
  theGsfELabel       = pset.getParameter<edm::InputTag>("ElectronCollectionLabel");

  // Minimum Pt for leptons for skimming
  stiffMinPt         = pset.getParameter<double>("stiffMinimumPt");
  softMinPt          = pset.getParameter<double>("softMinimumPt");
  nStiffLeptonMin    = pset.getParameter<int>("nStiffLeptonMinimum");
  nLeptonMin         = pset.getParameter<int>("nLeptonMinimum");

  nEvents         = 0;
  nSelectedEvents = 0;

}
HiggsToZZ4LeptonsSkim::~HiggsToZZ4LeptonsSkim ( )

Definition at line 57 of file HiggsToZZ4LeptonsSkim.cc.

References nEvents.

                                              {

  edm::LogVerbatim("HiggsToZZ4LeptonsSkim") 
  << " Number_events_read " << nEvents          
  << " Number_events_kept " << nSelectedEvents 
  << " Efficiency         " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl;
}

Member Function Documentation

bool HiggsToZZ4LeptonsSkim::filter ( edm::Event event,
const edm::EventSetup setup 
) [virtual]

Get event properties to send to builder to fill seed collection.

Implements edm::EDFilter.

Definition at line 68 of file HiggsToZZ4LeptonsSkim.cc.

References HI_PhotonSkim_cff::electrons, edm::HandleBase::isValid(), patZpeak::muons, nEvents, and edm::Handle< T >::product().

                                                                               {

  nEvents++;

  using reco::TrackCollection;

  bool keepEvent   = false;
  int  nStiffLeptons    = 0;
  int  nLeptons    = 0;
  

  // First look at muons:

  // Get the muon track collection from the event
  edm::Handle<reco::TrackCollection> muTracks;
  event.getByLabel(theGLBMuonLabel.label(), muTracks);

  if ( muTracks.isValid() ) {  
  
    reco::TrackCollection::const_iterator muons;
        
    // Loop over muon collections and count how many muons there are, 
    // and how many are above threshold
    for ( muons = muTracks->begin(); muons != muTracks->end(); ++muons ) {
      if ( muons->pt() > stiffMinPt) nStiffLeptons++; 
      if ( muons->pt() > softMinPt) nLeptons++; 
    }  
  } 
  
  // Now look at electrons:

  // Get the electron track collection from the event
  edm::Handle<reco::GsfElectronCollection> pTracks;

  event.getByLabel(theGsfELabel.label(),pTracks);

  if ( pTracks.isValid() ) {  

    const reco::GsfElectronCollection* eTracks = pTracks.product();

    reco::GsfElectronCollection::const_iterator electrons;

    // Loop over electron collections and count how many muons there are, 
    // and how many are above threshold
    for ( electrons = eTracks->begin(); electrons != eTracks->end(); ++electrons ) {
      float pt_e = electrons->pt(); 
      if ( pt_e > stiffMinPt) nStiffLeptons++; 
      if ( pt_e > softMinPt) nLeptons++; 
    }
  }

  // Make decision:
  if ( nStiffLeptons >= nStiffLeptonMin && nLeptons >= nLeptonMin) keepEvent = true;

  if (keepEvent) nSelectedEvents++;

  return keepEvent;
}

Member Data Documentation

Definition at line 45 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 42 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 49 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 42 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 48 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 52 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 47 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 46 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 53 of file HiggsToZZ4LeptonsSkim.h.

Definition at line 54 of file HiggsToZZ4LeptonsSkim.h.