#include <HiggsTo2GammaSkim.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
Get event properties to send to builder to fill seed collection. | |
HiggsTo2GammaSkim (const edm::ParameterSet &) | |
~HiggsTo2GammaSkim () | |
Private Attributes | |
bool | debug |
int | nEvents |
int | nPhotonMin |
int | nSelectedEvents |
float | photon1MinPt |
float | photon2MinPt |
edm::InputTag | thePhotonLabel |
Definition at line 26 of file HiggsTo2GammaSkim.h.
HiggsTo2GammaSkim::HiggsTo2GammaSkim | ( | const edm::ParameterSet & | pset | ) | [explicit] |
Definition at line 34 of file HiggsTo2GammaSkim.cc.
References debug, edm::ParameterSet::getParameter(), and nEvents.
{ // Local Debug flag debug = pset.getParameter<bool>("DebugHiggsTo2GammaSkim"); // Reconstructed objects thePhotonLabel = pset.getParameter<edm::InputTag>("PhotonCollectionLabel"); // Minimum Pt for photons for skimming photon1MinPt = pset.getParameter<double>("photon1MinimumPt"); nPhotonMin = pset.getParameter<int>("nPhotonMinimum"); nEvents = 0; nSelectedEvents = 0; }
HiggsTo2GammaSkim::~HiggsTo2GammaSkim | ( | ) |
Definition at line 54 of file HiggsTo2GammaSkim.cc.
References nEvents.
{ edm::LogVerbatim("HiggsTo2GammaSkim") << " Number_events_read " << nEvents << " Number_events_kept " << nSelectedEvents << " Efficiency " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl; }
bool HiggsTo2GammaSkim::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 65 of file HiggsTo2GammaSkim.cc.
References edm::HandleBase::isValid(), nEvents, interactiveExample::photons, and edm::Handle< T >::product().
{ nEvents++; using reco::PhotonCollection; bool keepEvent = false; int nPhotons = 0; // Look at photons: // Get the photon collection from the event edm::Handle<reco::PhotonCollection> photonHandle; event.getByLabel(thePhotonLabel.label(),photonHandle); if ( photonHandle.isValid() ) { const reco::PhotonCollection* phoCollection = photonHandle.product(); reco::PhotonCollection::const_iterator photons; // Loop over photon collections and count how many photons there are, // and how many are above the thresholds // Question: do we need to take the reconstructed primary vertex at this point? // Here, I assume that the et is taken with respect to the nominal vertex (0,0,0). for ( photons = phoCollection->begin(); photons != phoCollection->end(); ++photons ) { float et_p = photons->et(); if ( et_p > photon1MinPt) nPhotons++; } } // Make decision: if ( nPhotons >= nPhotonMin ) keepEvent = true; if (keepEvent) nSelectedEvents++; return keepEvent; }
bool HiggsTo2GammaSkim::debug [private] |
Definition at line 43 of file HiggsTo2GammaSkim.h.
int HiggsTo2GammaSkim::nEvents [private] |
Definition at line 40 of file HiggsTo2GammaSkim.h.
int HiggsTo2GammaSkim::nPhotonMin [private] |
Definition at line 46 of file HiggsTo2GammaSkim.h.
int HiggsTo2GammaSkim::nSelectedEvents [private] |
Definition at line 40 of file HiggsTo2GammaSkim.h.
float HiggsTo2GammaSkim::photon1MinPt [private] |
Definition at line 44 of file HiggsTo2GammaSkim.h.
float HiggsTo2GammaSkim::photon2MinPt [private] |
Definition at line 45 of file HiggsTo2GammaSkim.h.
Definition at line 49 of file HiggsTo2GammaSkim.h.