#include <HLTSmartSinglet.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HLTSmartSinglet (const edm::ParameterSet &) | |
~HLTSmartSinglet () | |
Private Attributes | |
std::string | cut_ |
edm::InputTag | inputTag_ |
int | min_N_ |
bool | saveTags_ |
StringCutObjectSelector< T, true > | select_ |
This class is an HLTFilter (-> EDFilter) implementing a smart HLT trigger cut, specified as a string such as "pt>15 && -3<eta<3", for single objects of the same physics type, allowing to cut on variables relating to both the base class T and the derived actual class
See header file for documentation
Definition at line 31 of file HLTSmartSinglet.h.
HLTSmartSinglet< T, Tid >::HLTSmartSinglet | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 26 of file HLTSmartSinglet.cc.
References HLTSmartSinglet< T, Tid >::cut_, edm::InputTag::encode(), HLTSmartSinglet< T, Tid >::inputTag_, LogDebug, and HLTSmartSinglet< T, Tid >::min_N_.
: inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")), saveTags_ (iConfig.template getParameter<bool>("saveTags")), cut_ (iConfig.template getParameter<std::string> ("cut" )), min_N_ (iConfig.template getParameter<int> ("MinN" )), select_ (cut_ ) { LogDebug("") << "Input/cut/ncut : " << inputTag_.encode() << " " << cut_<< " " << min_N_ ; //register your products produces<trigger::TriggerFilterObjectWithRefs>(); }
HLTSmartSinglet< T, Tid >::~HLTSmartSinglet | ( | ) |
Definition at line 40 of file HLTSmartSinglet.cc.
{ }
bool HLTSmartSinglet< T, Tid >::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements HLTFilter.
Definition at line 51 of file HLTSmartSinglet.cc.
References accept(), edm::Event::getByLabel(), i, module(), n, path(), edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, and MultipleCompare::ref.
{ 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 (select_(*i)) { n++; ref=TRef(objects,distance(objects->begin(),i)); filterObject->addObject(Tid,ref); } } // filter decision bool accept(n>=min_N_); // put filter object into the Event iEvent.put(filterObject); return accept; }
std::string HLTSmartSinglet< T, Tid >::cut_ [private] |
Definition at line 42 of file HLTSmartSinglet.h.
Referenced by HLTSmartSinglet< T, Tid >::HLTSmartSinglet().
edm::InputTag HLTSmartSinglet< T, Tid >::inputTag_ [private] |
Definition at line 40 of file HLTSmartSinglet.h.
Referenced by HLTSmartSinglet< T, Tid >::HLTSmartSinglet().
int HLTSmartSinglet< T, Tid >::min_N_ [private] |
Definition at line 43 of file HLTSmartSinglet.h.
Referenced by HLTSmartSinglet< T, Tid >::HLTSmartSinglet().
bool HLTSmartSinglet< T, Tid >::saveTags_ [private] |
Definition at line 41 of file HLTSmartSinglet.h.
StringCutObjectSelector<T,true> HLTSmartSinglet< T, Tid >::select_ [private] |
Definition at line 45 of file HLTSmartSinglet.h.