#include <HLTSinglet.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HLTSinglet (const edm::ParameterSet &) | |
~HLTSinglet () | |
Private Attributes | |
edm::InputTag | inputTag_ |
double | max_Eta_ |
int | min_N_ |
double | min_Pt_ |
bool | saveTags_ |
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 27 of file HLTSinglet.h.
HLTSinglet< T, Tid >::HLTSinglet | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 77 of file HLTSinglet.cc.
References edm::InputTag::encode(), HLTSinglet< T, Tid >::inputTag_, LogDebug, HLTSinglet< T, Tid >::max_Eta_, HLTSinglet< T, Tid >::min_N_, and HLTSinglet< T, Tid >::min_Pt_.
: inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")), saveTags_ (iConfig.template getParameter<bool>("saveTags")), min_Pt_ (iConfig.template getParameter<double> ("MinPt" )), max_Eta_ (iConfig.template getParameter<double> ("MaxEta" )), min_N_ (iConfig.template getParameter<int> ("MinN" )) { LogDebug("") << "Input/ptcut/etacut/ncut : " << inputTag_.encode() << " " << min_Pt_ << " " << max_Eta_ << " " << min_N_ ; //register your products produces<trigger::TriggerFilterObjectWithRefs>(); }
HLTSinglet< T, Tid >::~HLTSinglet | ( | ) |
Definition at line 93 of file HLTSinglet.cc.
{ }
bool HLTSinglet< T, Tid >::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements HLTFilter.
Definition at line 104 of file HLTSinglet.cc.
References abs, accept(), edm::Event::getByLabel(), i, module(), n, path(), edm::Event::put(), and dt_dqm_sourceclient_common_cff::reco.
{ 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 auto_ptr<TriggerFilterObjectWithRefs> filterobject (new TriggerFilterObjectWithRefs(path(),module())); if (saveTags_) filterobject->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->pt() >= min_Pt_) && ( (max_Eta_ < 0.0) || (std::abs(i->eta()) <= max_Eta_) ) ) { n++; ref=TRef(objects,distance(objects->begin(),i)); filterobject->addObject(getObjectType<T, Tid>(*i),ref); } } // filter decision bool accept(n>=min_N_); // put filter object into the Event iEvent.put(filterobject); return accept; }
edm::InputTag HLTSinglet< T, Tid >::inputTag_ [private] |
Definition at line 35 of file HLTSinglet.h.
Referenced by HLTSinglet< T, Tid >::HLTSinglet().
double HLTSinglet< T, Tid >::max_Eta_ [private] |
Definition at line 38 of file HLTSinglet.h.
Referenced by HLTSinglet< T, Tid >::HLTSinglet().
int HLTSinglet< T, Tid >::min_N_ [private] |
Definition at line 39 of file HLTSinglet.h.
Referenced by HLTSinglet< T, Tid >::HLTSinglet().
double HLTSinglet< T, Tid >::min_Pt_ [private] |
Definition at line 37 of file HLTSinglet.h.
Referenced by HLTSinglet< T, Tid >::HLTSinglet().
bool HLTSinglet< T, Tid >::saveTags_ [private] |
Definition at line 36 of file HLTSinglet.h.