![]() |
![]() |
#include <HLTForwardBackwardJetsFilter.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HLTForwardBackwardJetsFilter (const edm::ParameterSet &) | |
~HLTForwardBackwardJetsFilter () | |
Static Public Member Functions | |
static void | fillDescriptions (edm::ConfigurationDescriptions &descriptions) |
Private Attributes | |
edm::InputTag | inputTag_ |
double | maxEta_ |
double | minEta_ |
double | minPt_ |
bool | saveTag_ |
Definition at line 18 of file HLTForwardBackwardJetsFilter.h.
HLTForwardBackwardJetsFilter::HLTForwardBackwardJetsFilter | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 28 of file HLTForwardBackwardJetsFilter.cc.
References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), inputTag_, maxEta_, minEta_, minPt_, and saveTag_.
{ inputTag_ = iConfig.getParameter< edm::InputTag > ("inputTag"); saveTag_ = iConfig.getUntrackedParameter<bool>("saveTag"); minPt_ = iConfig.getParameter<double> ("minPt"); minEta_ = iConfig.getParameter<double> ("minEta"); maxEta_ = iConfig.getParameter<double> ("maxEta"); //register your products produces<trigger::TriggerFilterObjectWithRefs>(); }
HLTForwardBackwardJetsFilter::~HLTForwardBackwardJetsFilter | ( | ) |
Definition at line 40 of file HLTForwardBackwardJetsFilter.cc.
{}
void HLTForwardBackwardJetsFilter::fillDescriptions | ( | edm::ConfigurationDescriptions & | descriptions | ) | [static] |
Reimplemented from edm::EDFilter.
Definition at line 43 of file HLTForwardBackwardJetsFilter.cc.
References edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::add(), and edm::ParameterSetDescription::addUntracked().
{ edm::ParameterSetDescription desc; desc.add<edm::InputTag>("inputTag",edm::InputTag("hltIterativeCone5CaloJetsRegional")); desc.addUntracked<bool>("saveTag",false); desc.add<double>("minPt",15.0); desc.add<double>("minEta",3.0); desc.add<double>("maxEta",5.1); descriptions.add("hltForwardBackwardJetsFilter",desc); }
bool HLTForwardBackwardJetsFilter::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements HLTFilter.
Definition at line 55 of file HLTForwardBackwardJetsFilter.cc.
References accept(), edm::Event::getByLabel(), inputTag_, maxEta_, minEta_, minPt_, module(), path(), edm::Event::put(), saveTag_, and trigger::TriggerJet.
{ using namespace trigger; // The filter object std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterobject (new trigger::TriggerFilterObjectWithRefs(path(),module())); if (saveTag_) filterobject->addCollectionTag(inputTag_); edm::Handle<reco::CaloJetCollection> recocalojets; iEvent.getByLabel(inputTag_,recocalojets); // look at all candidates, check cuts and add to filter object int nplusjets(0); int nminusjets(0); if(recocalojets->size() > 1){ // events with two or more jets // look for jets satifying pt and eta cuts; first on the plus side, then the minus side for (reco::CaloJetCollection::const_iterator recocalojet = recocalojets->begin(); recocalojet!=(recocalojets->end()); recocalojet++) { float ptjet=recocalojet->pt(); float etajet=recocalojet->eta(); if( ptjet > minPt_ ){ if ( etajet > minEta_ && etajet < maxEta_ ){ nplusjets++; reco::CaloJetRef ref(reco::CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet))); filterobject->addObject(TriggerJet,ref); } } } if (nplusjets > 0) { for (reco::CaloJetCollection::const_iterator recocalojet = recocalojets->begin(); recocalojet!=(recocalojets->end()); recocalojet++) { float ptjet=recocalojet->pt(); float etajet=recocalojet->eta(); if( ptjet > minPt_ ){ if ( etajet < -minEta_ && etajet > -maxEta_ ){ nminusjets++; reco::CaloJetRef ref(reco::CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet))); filterobject->addObject(TriggerJet,ref); } } } } } // events with two or more jets // filter decision bool accept(nplusjets>0 && nminusjets>0); // put filter object into the Event iEvent.put(filterobject); return accept; }
Definition at line 27 of file HLTForwardBackwardJetsFilter.h.
Referenced by filter(), and HLTForwardBackwardJetsFilter().
double HLTForwardBackwardJetsFilter::maxEta_ [private] |
Definition at line 31 of file HLTForwardBackwardJetsFilter.h.
Referenced by filter(), and HLTForwardBackwardJetsFilter().
double HLTForwardBackwardJetsFilter::minEta_ [private] |
Definition at line 30 of file HLTForwardBackwardJetsFilter.h.
Referenced by filter(), and HLTForwardBackwardJetsFilter().
double HLTForwardBackwardJetsFilter::minPt_ [private] |
Definition at line 29 of file HLTForwardBackwardJetsFilter.h.
Referenced by filter(), and HLTForwardBackwardJetsFilter().
bool HLTForwardBackwardJetsFilter::saveTag_ [private] |
Definition at line 28 of file HLTForwardBackwardJetsFilter.h.
Referenced by filter(), and HLTForwardBackwardJetsFilter().