CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTForwardBackwardJetsFilter.cc
Go to the documentation of this file.
1 
9 
13 
21 
22 #include<typeinfo>
23 
24 //
25 // constructors and destructor
26 //
27 template<typename T>
29  HLTFilter(iConfig),
30  inputTag_ (iConfig.template getParameter< edm::InputTag > ("inputTag")),
31  minPt_ (iConfig.template getParameter<double> ("minPt")),
32  minEta_ (iConfig.template getParameter<double> ("minEta")),
33  maxEta_ (iConfig.template getParameter<double> ("maxEta")),
34  nNeg_ (iConfig.template getParameter<unsigned int>("nNeg")),
35  nPos_ (iConfig.template getParameter<unsigned int>("nPos")),
36  nTot_ (iConfig.template getParameter<unsigned int>("nTot")),
37  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
38 {
39  LogDebug("") << "HLTForwardBackwardJetsFilter: Input/minPt/minEta/maxEta/triggerType : "
40  << inputTag_.encode() << " "
41  << minPt_ << " "
42  << minEta_ << " "
43  << maxEta_ << " "
44  << nNeg_ << " "
45  << nPos_ << " "
46  << nTot_ << " "
47  << triggerType_;
48 }
49 
50 template<typename T>
52 
53 template<typename T>
54 void
57  makeHLTFilterDescription(desc);
58  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltIterativeCone5CaloJetsRegional"));
59  desc.add<double>("minPt",15.0);
60  desc.add<double>("minEta",3.0);
61  desc.add<double>("maxEta",5.1);
62  desc.add<unsigned int>("nNeg",1);
63  desc.add<unsigned int>("nPos",1);
64  desc.add<unsigned int>("nTot",0);
65  desc.add<int>("triggerType",trigger::TriggerJet);
66  descriptions.add(std::string("hlt")+std::string(typeid(HLTForwardBackwardJetsFilter<T>).name()),desc);
67 }
68 
69 // ------------ method called to produce the data ------------
70 template<typename T>
71 bool
73 {
74  using namespace std;
75  using namespace edm;
76  using namespace reco;
77  using namespace trigger;
78 
79  typedef vector<T> TCollection;
80  typedef Ref<TCollection> TRef;
81 
82  // The filter object
83  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
84 
85  // get hold of collection of objects
86  Handle<TCollection> objects;
87  iEvent.getByLabel(inputTag_,objects);
88 
89  // look at all candidates, check cuts and add to filter object
90  unsigned int nPosJets(0);
91  unsigned int nNegJets(0);
92 
93  typename TCollection::const_iterator jet;
94  // look for jets satifying pt and eta cuts; first on the plus side, then the minus side
95 
96  for (jet=objects->begin(); jet!=objects->end(); jet++) {
97  double ptjet = jet->pt();
98  double etajet = jet->eta();
99  if( ptjet >= minPt_ ){
100  if (( minEta_<= etajet) && (etajet <= maxEta_) ){
101  nPosJets++;
102  TRef ref = TRef(objects,distance(objects->begin(),jet));
103  filterproduct.addObject(triggerType_,ref);
104  }
105  if ((-maxEta_<= etajet) && (etajet <=-minEta_) ){
106  nNegJets++;
107  TRef ref = TRef(objects,distance(objects->begin(),jet));
108  filterproduct.addObject(triggerType_,ref);
109  }
110  }
111  }
112 
113  // filter decision
114  const bool accept(
115  ( nNegJets >= nNeg_ ) &&
116  ( nPosJets >= nPos_ ) &&
117  ((nNegJets+nPosJets) >= nTot_ )
118  );
119 
120  return accept;
121 }
#define LogDebug(id)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
std::string encode() const
Definition: InputTag.cc:72
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
int iEvent
Definition: GenABIO.cc:243
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLTForwardBackwardJetsFilter(const edm::ParameterSet &)
def template
Definition: svgfig.py:520