#include <HLTGlobalSums.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HLTGlobalSums (const edm::ParameterSet &) | |
~HLTGlobalSums () | |
Private Attributes | |
edm::InputTag | inputTag_ |
double | max_ |
double | min_ |
int | min_N_ |
std::string | observable_ |
bool | saveTag_ |
int | tid_ |
This class is an HLTFilter (-> EDFilter) implementing cuts on global sums such as the scalar sum of Et (a.k.a. H_T), available in the T=CaloMET or T=MET object.
Definition at line 26 of file HLTGlobalSums.h.
HLTGlobalSums< T, Tid >::HLTGlobalSums | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 28 of file HLTGlobalSums.cc.
References edm::InputTag::encode(), HLTGlobalSums< T, Tid >::inputTag_, LogDebug, HLTGlobalSums< T, Tid >::max_, HLTGlobalSums< T, Tid >::min_, HLTGlobalSums< T, Tid >::min_N_, HLTGlobalSums< T, Tid >::observable_, HLTGlobalSums< T, Tid >::tid_, trigger::TriggerELongit, trigger::TriggerHLongit, trigger::TriggerMETSig, trigger::TriggerMHTSig, trigger::TriggerTET, and trigger::TriggerTHT.
: inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")), saveTag_ (iConfig.template getUntrackedParameter<bool> ("saveTag",false)), observable_ (iConfig.template getParameter<std::string>("observable")), min_ (iConfig.template getParameter<double>("Min")), max_ (iConfig.template getParameter<double>("Max")), min_N_ (iConfig.template getParameter<int>("MinN")), tid_() { LogDebug("") << "InputTags and cuts : " << inputTag_.encode() << " " << observable_ << " Range [" << min_ << " " << max_ << "]" << " MinN =" << min_N_ ; if (observable_=="sumEt") { tid_=Tid; } else if (observable_=="mEtSig") { if (Tid==trigger::TriggerTET) { tid_=trigger::TriggerMETSig; } else if (Tid==trigger::TriggerTHT) { tid_=trigger::TriggerMHTSig; } else { tid_=Tid; } } else if (observable_=="e_longitudinal") { if (Tid==trigger::TriggerTET) { tid_=trigger::TriggerELongit; } else if (Tid==trigger::TriggerTHT) { tid_=trigger::TriggerHLongit; } else { tid_=Tid; } } else { tid_=Tid; } //register your products produces<trigger::TriggerFilterObjectWithRefs>(); }
HLTGlobalSums< T, Tid >::~HLTGlobalSums | ( | ) |
Definition at line 70 of file HLTGlobalSums.cc.
{ }
bool HLTGlobalSums< T, Tid >::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements HLTFilter.
Definition at line 81 of file HLTGlobalSums.cc.
References abs, accept(), edm::Event::getByLabel(), LogDebug, module(), n, path(), edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, trigger::TriggerELongit, trigger::TriggerHLongit, trigger::TriggerMETSig, trigger::TriggerMHTSig, trigger::TriggerTET, trigger::TriggerTHT, and relativeConstraints::value.
{ 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 (saveTag_) filterobject->addCollectionTag(inputTag_); // Ref to Candidate object to be recorded in filter object TRef ref; // get hold of MET product from Event Handle<TCollection> objects; iEvent.getByLabel(inputTag_,objects); if (!objects.isValid()) { LogDebug("") << inputTag_ << " collection not found!"; iEvent.put(filterobject); return false; } LogDebug("") << "Size of MET collection: " << objects->size(); if (objects->size()==0) { LogDebug("") << "MET collection does not contain a MET object!"; } else if (objects->size()>1) { LogDebug("") << "MET collection contains more than one MET object!"; } int n(0); double value(0.0); typename TCollection::const_iterator ibegin(objects->begin()); typename TCollection::const_iterator iend(objects->end()); typename TCollection::const_iterator iter; for (iter=ibegin; iter!=iend; iter++) { // get hold of value of observable to cut on if ( (tid_==TriggerTET) || (tid_==TriggerTHT) ) { value=iter->sumEt(); } else if ( (tid_==TriggerMETSig) || (tid_==TriggerMHTSig) ) { value=iter->mEtSig(); } else if ( (tid_==TriggerELongit) || (tid_==TriggerHLongit) ) { value=iter->e_longitudinal(); } else { value=0.0; } value=std::abs(value); if ( ( (min_<0.0) || (min_<=value) ) && ( (max_<0.0) || (value<=max_) ) ) { n++; ref=TRef(objects,distance(ibegin,iter)); filterobject->addObject(tid_,ref); } } // filter decision const bool accept(n>=min_N_); // put filter object into the Event iEvent.put(filterobject); return accept; }
edm::InputTag HLTGlobalSums< T, Tid >::inputTag_ [private] |
Definition at line 36 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::HLTGlobalSums().
double HLTGlobalSums< T, Tid >::max_ [private] |
Definition at line 39 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::HLTGlobalSums().
double HLTGlobalSums< T, Tid >::min_ [private] |
Definition at line 39 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::HLTGlobalSums().
int HLTGlobalSums< T, Tid >::min_N_ [private] |
Definition at line 40 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::HLTGlobalSums().
std::string HLTGlobalSums< T, Tid >::observable_ [private] |
Definition at line 38 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::HLTGlobalSums().
bool HLTGlobalSums< T, Tid >::saveTag_ [private] |
Definition at line 37 of file HLTGlobalSums.h.
int HLTGlobalSums< T, Tid >::tid_ [private] |
Definition at line 41 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::HLTGlobalSums().