CMS 3D CMS Logo

HLTGlobalSums.cc
Go to the documentation of this file.
1 
10 #include <cmath>
11 
13 #include "HLTGlobalSums.h"
14 
16 
19 
21 
23 
24 //
25 // constructors and destructor
26 //
27 template<typename T>
29  inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")),
30  inputToken_ (consumes<std::vector<T> >(inputTag_)),
31  triggerType_(iConfig.template getParameter<int>("triggerType")),
32  observable_ (iConfig.template getParameter<std::string>("observable")),
33  min_ (iConfig.template getParameter<double>("Min")),
34  max_ (iConfig.template getParameter<double>("Max")),
35  min_N_ (iConfig.template getParameter<int>("MinN")),
36  tid_(triggerType_)
37 {
38  LogDebug("") << "InputTags and cuts : "
39  << inputTag_.encode() << " "
40  << triggerType_ << " "
41  << observable_
42  << " Range [" << min_ << " " << max_ << "]"
43  << " MinN =" << min_N_ ;
44 
45  if (observable_=="sumEt") {
47  } else if (observable_=="mEtSig") {
50  } else if (triggerType_==trigger::TriggerTHT) {
52  } else {
54  }
55  } else if (observable_=="e_longitudinal") {
58  } else if (triggerType_==trigger::TriggerTHT) {
60  } else {
62  }
63  } else {
65  }
66 }
67 
68 template<typename T>
70 
71 template<typename T>
72 void
76  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltCollection"));
77  desc.add<int>("triggerType",0);
78  desc.add<std::string>("observable","");
79  desc.add<double>("Min",-1e125);
80  desc.add<double>("Max",+1e125);
81  desc.add<int>("MinN",1);
82  descriptions.add(defaultModuleLabel<HLTGlobalSums<T>>(), desc);
83 }
84 
85 //
86 // member functions
87 //
88 
89 // ------------ method called to produce the data ------------
90 template<typename T>
91 bool
93 {
94  using namespace std;
95  using namespace edm;
96  using namespace reco;
97  using namespace trigger;
98 
99  typedef vector<T> TCollection;
100  typedef Ref<TCollection> TRef;
101 
102  // All HLT filters must create and fill an HLT filter object,
103  // recording any reconstructed physics objects satisfying (or not)
104  // this HLT filter, and place it in the Event.
105 
106  // The filter object
107  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
108  // Ref to Candidate object to be recorded in filter object
109  TRef ref;
110 
111 
112  // get hold of MET product from Event
114  iEvent.getByToken(inputToken_,objects);
115  if (!objects.isValid()) {
116  LogDebug("") << inputTag_ << " collection not found!";
117  return false;
118  }
119 
120  LogDebug("") << "Size of MET collection: " << objects->size();
121  if (objects->empty()) {
122  LogDebug("") << "MET collection does not contain a MET object!";
123  } else if (objects->size()>1) {
124  LogDebug("") << "MET collection contains more than one MET object!";
125  }
126 
127  int n(0);
128  double value(0.0);
129  typename TCollection::const_iterator ibegin(objects->begin());
130  typename TCollection::const_iterator iend(objects->end());
131  typename TCollection::const_iterator iter;
132  for (iter=ibegin; iter!=iend; iter++) {
133 
134  // get hold of value of observable to cut on
135  if ( (tid_==TriggerTET) || (tid_==TriggerTHT) ) {
136  value=iter->sumEt();
137  } else if ( (tid_==TriggerMETSig) || (tid_==TriggerMHTSig) ) {
138  value=iter->mEtSig();
139  } else if ( (tid_==TriggerELongit) || (tid_==TriggerHLongit) ) {
140  value=iter->e_longitudinal();
141  } else {
142  value=0.0;
143  }
144 
145  value=std::abs(value);
146 
147  if ( ( (min_<0.0) || (min_<=value) ) &&
148  ( (max_<0.0) || (value<=max_) ) ) {
149  n++;
150  ref=TRef(objects,distance(ibegin,iter));
151  filterproduct.addObject(tid_,ref);
152  }
153 
154  }
155 
156  // filter decision
157  const bool accept(n>=min_N_);
158 
159  return accept;
160 }
#define LogDebug(id)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
edm::EDGetTokenT< std::vector< T > > inputToken_
Definition: HLTGlobalSums.h:38
std::string defaultModuleLabel()
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::string encode() const
Definition: InputTag.cc:166
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:520
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
~HLTGlobalSums() override
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
edm::InputTag inputTag_
Definition: HLTGlobalSums.h:37
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
long double T
HLTGlobalSums(const edm::ParameterSet &)
std::string observable_
Definition: HLTGlobalSums.h:40
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override