#include <HLTSinglet.h>
Public Member Functions | |
virtual bool | hltFilter (edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) |
HLTSinglet (const edm::ParameterSet &) | |
~HLTSinglet () | |
Static Public Member Functions | |
static void | fillDescriptions (edm::ConfigurationDescriptions &descriptions) |
Private Attributes | |
edm::InputTag | inputTag_ |
double | max_Eta_ |
double | min_E_ |
double | min_Mass_ |
int | min_N_ |
double | min_Pt_ |
int | tid_ |
int | triggerType_ |
This class is an HLTFilter (-> EDFilter) implementing a basic HLT trigger for single objects of the same physics type, cutting on variables relating to their 4-momentum representation
See header file for documentation
Definition at line 28 of file HLTSinglet.h.
HLTSinglet< T >::HLTSinglet | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 78 of file HLTSinglet.cc.
References edm::InputTag::encode(), HLTSinglet< T >::inputTag_, LogDebug, HLTSinglet< T >::max_Eta_, HLTSinglet< T >::min_E_, HLTSinglet< T >::min_Mass_, HLTSinglet< T >::min_N_, and HLTSinglet< T >::min_Pt_.
: HLTFilter(iConfig), inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")), triggerType_ (iConfig.template getParameter<int>("triggerType")), min_E_ (iConfig.template getParameter<double> ("MinE" )), min_Pt_ (iConfig.template getParameter<double> ("MinPt" )), min_Mass_ (iConfig.template getParameter<double> ("MinMass" )), max_Eta_ (iConfig.template getParameter<double> ("MaxEta" )), min_N_ (iConfig.template getParameter<int> ("MinN" )), tid_ (triggerType_) { LogDebug("") << "Input/ptcut/etacut/ncut : " << inputTag_.encode() << " " << min_E_ << " " << min_Pt_ << " " << min_Mass_ << " " << max_Eta_ << " " << min_N_ ; }
HLTSinglet< T >::~HLTSinglet | ( | ) |
Definition at line 95 of file HLTSinglet.cc.
{ }
void HLTSinglet< T >::fillDescriptions | ( | edm::ConfigurationDescriptions & | descriptions | ) | [static] |
Reimplemented from edm::EDFilter.
Definition at line 101 of file HLTSinglet.cc.
References edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::add(), and mergeVDriftHistosByStation::name.
{ edm::ParameterSetDescription desc; makeHLTFilterDescription(desc); desc.add<edm::InputTag>("inputTag",edm::InputTag("hltCollection")); desc.add<int>("triggerType",0); desc.add<double>("MinE",-1.0); desc.add<double>("MinPt",-1.0); desc.add<double>("MinMass",-1.0); desc.add<double>("MaxEta",-1.0); desc.add<int>("MinN",1); descriptions.add(std::string("hlt")+std::string(typeid(HLTSinglet<T>).name()),desc); }
bool HLTSinglet< T >::hltFilter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup, | ||
trigger::TriggerFilterObjectWithRefs & | filterproduct | ||
) | [virtual] |
Implements HLTFilter.
Definition at line 121 of file HLTSinglet.cc.
References abs, accept(), trigger::TriggerFilterObjectWithRefs::addCollectionTag(), trigger::TriggerRefsCollections::addObject(), edm::Event::getByLabel(), i, n, dt_dqm_sourceclient_common_cff::reco, and CommPDSkim_cfg::saveTags.
{ using namespace std; using namespace edm; using namespace reco; using namespace trigger; typedef vector<T> TCollection; typedef Ref<TCollection> TRef; // All HLT filters must create and fill an HLT filter object, // recording any reconstructed physics objects satisfying (or not) // this HLT filter, and place it in the Event. // The filter object if (saveTags()) filterproduct.addCollectionTag(inputTag_); // Ref to Candidate object to be recorded in filter object TRef ref; // get hold of collection of objects Handle<TCollection> objects; iEvent.getByLabel (inputTag_,objects); // look at all objects, check cuts and add to filter object int n(0); typename TCollection::const_iterator i ( objects->begin() ); for (; i!=objects->end(); i++) { if ( (i->energy() >= min_E_) && (i->pt() >= min_Pt_) && (i->mass() >= min_Mass_) && ( (max_Eta_ < 0.0) || (std::abs(i->eta()) <= max_Eta_) ) ) { n++; ref=TRef(objects,distance(objects->begin(),i)); tid_=getObjectType<T>(*i); if (tid_==0) tid_=triggerType_; filterproduct.addObject(tid_,ref); } } // filter decision bool accept(n>=min_N_); return accept; }
edm::InputTag HLTSinglet< T >::inputTag_ [private] |
Definition at line 37 of file HLTSinglet.h.
Referenced by HLTSinglet< T >::HLTSinglet().
double HLTSinglet< T >::max_Eta_ [private] |
Definition at line 42 of file HLTSinglet.h.
Referenced by HLTSinglet< T >::HLTSinglet().
double HLTSinglet< T >::min_E_ [private] |
Definition at line 39 of file HLTSinglet.h.
Referenced by HLTSinglet< T >::HLTSinglet().
double HLTSinglet< T >::min_Mass_ [private] |
Definition at line 41 of file HLTSinglet.h.
Referenced by HLTSinglet< T >::HLTSinglet().
int HLTSinglet< T >::min_N_ [private] |
Definition at line 43 of file HLTSinglet.h.
Referenced by HLTSinglet< T >::HLTSinglet().
double HLTSinglet< T >::min_Pt_ [private] |
Definition at line 40 of file HLTSinglet.h.
Referenced by HLTSinglet< T >::HLTSinglet().
int HLTSinglet< T >::tid_ [private] |
Definition at line 44 of file HLTSinglet.h.
int HLTSinglet< T >::triggerType_ [private] |
Definition at line 38 of file HLTSinglet.h.