#include <HLTrigger/HLTfilters/interface/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_ |
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 | ) | [inline, 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::TriggerHT, and trigger::TriggerMETSig.
00028 : 00029 inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")), 00030 saveTag_ (iConfig.template getUntrackedParameter<bool> ("saveTag",false)), 00031 observable_ (iConfig.template getParameter<std::string>("observable")), 00032 min_ (iConfig.template getParameter<double>("Min")), 00033 max_ (iConfig.template getParameter<double>("Max")), 00034 min_N_ (iConfig.template getParameter<int>("MinN")), 00035 tid_() 00036 { 00037 LogDebug("") << "InputTags and cuts : " 00038 << inputTag_.encode() << " " << observable_ 00039 << " Range [" << min_ << " " << max_ << "]" 00040 << " MinN =" << min_N_ 00041 ; 00042 00043 if (observable_=="sumEt") { 00044 tid_=trigger::TriggerHT; 00045 } else if (observable_=="mEtSig") { 00046 tid_=trigger::TriggerMETSig; 00047 } else if (observable_=="e_longitudinal") { 00048 tid_=trigger::TriggerELongit; 00049 } else { 00050 tid_=Tid; 00051 } 00052 00053 //register your products 00054 produces<trigger::TriggerFilterObjectWithRefs>(); 00055 }
HLTGlobalSums< T, Tid >::~HLTGlobalSums | ( | ) | [inline] |
bool HLTGlobalSums< T, Tid >::filter | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [inline, virtual] |
Implements HLTFilter.
Definition at line 69 of file HLTGlobalSums.cc.
References funct::abs(), edm::Event::getByLabel(), HLTGlobalSums< T, Tid >::inputTag_, iter, LogDebug, HLTGlobalSums< T, Tid >::max_, HLTGlobalSums< T, Tid >::min_, HLTGlobalSums< T, Tid >::min_N_, module(), n, path(), edm::Event::put(), HcalSimpleRecAlgoImpl::reco(), HLTGlobalSums< T, Tid >::saveTag_, std, HLTGlobalSums< T, Tid >::tid_, trigger::TriggerELongit, trigger::TriggerHT, trigger::TriggerMETSig, and value.
00070 { 00071 using namespace std; 00072 using namespace edm; 00073 using namespace reco; 00074 using namespace trigger; 00075 00076 typedef vector<T> TCollection; 00077 typedef Ref<TCollection> TRef; 00078 00079 // All HLT filters must create and fill an HLT filter object, 00080 // recording any reconstructed physics objects satisfying (or not) 00081 // this HLT filter, and place it in the Event. 00082 00083 // The filter object 00084 auto_ptr<TriggerFilterObjectWithRefs> 00085 filterobject (new TriggerFilterObjectWithRefs(path(),module())); 00086 if (saveTag_) filterobject->addCollectionTag(inputTag_); 00087 // Ref to Candidate object to be recorded in filter object 00088 TRef ref; 00089 00090 00091 // get hold of MET product from Event 00092 Handle<TCollection> objects; 00093 iEvent.getByLabel(inputTag_,objects); 00094 if (!objects.isValid()) { 00095 LogDebug("") << inputTag_ << " collection not found!"; 00096 iEvent.put(filterobject); 00097 return false; 00098 } 00099 00100 LogDebug("") << "Size of MET collection: " << objects->size(); 00101 if (objects->size()==0) { 00102 LogDebug("") << "MET collection does not contain a MET object!"; 00103 } else if (objects->size()>1) { 00104 LogDebug("") << "MET collection contains more than one MET object!"; 00105 } 00106 00107 int n(0); 00108 double value(0.0); 00109 typename TCollection::const_iterator ibegin(objects->begin()); 00110 typename TCollection::const_iterator iend(objects->end()); 00111 typename TCollection::const_iterator iter; 00112 for (iter=ibegin; iter!=iend; iter++) { 00113 00114 // get hold of value of observable to cut on 00115 if (tid_==TriggerHT) { 00116 value=iter->sumEt(); 00117 } else if (tid_==TriggerMETSig) { 00118 value=iter->mEtSig(); 00119 } else if (tid_==TriggerELongit) { 00120 value=iter->e_longitudinal(); 00121 } else { 00122 value=0.0; 00123 } 00124 00125 value=abs(value); 00126 00127 if ( ( (min_<0.0) || (min_<=value) ) && 00128 ( (max_<0.0) || (value<=max_) ) ) { 00129 n++; 00130 ref=TRef(objects,distance(ibegin,iter)); 00131 filterobject->addObject(tid_,ref); 00132 } 00133 00134 } 00135 00136 // filter decision 00137 const bool accept(n>=min_N_); 00138 00139 // put filter object into the Event 00140 iEvent.put(filterobject); 00141 00142 return accept; 00143 }
edm::InputTag HLTGlobalSums< T, Tid >::inputTag_ [private] |
Definition at line 36 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::filter(), and HLTGlobalSums< T, Tid >::HLTGlobalSums().
double HLTGlobalSums< T, Tid >::max_ [private] |
Definition at line 39 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::filter(), and HLTGlobalSums< T, Tid >::HLTGlobalSums().
double HLTGlobalSums< T, Tid >::min_ [private] |
Definition at line 39 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::filter(), and HLTGlobalSums< T, Tid >::HLTGlobalSums().
int HLTGlobalSums< T, Tid >::min_N_ [private] |
Definition at line 40 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::filter(), and 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] |
int HLTGlobalSums< T, Tid >::tid_ [private] |
Definition at line 41 of file HLTGlobalSums.h.
Referenced by HLTGlobalSums< T, Tid >::filter(), and HLTGlobalSums< T, Tid >::HLTGlobalSums().