CMS 3D CMS Logo

HLT2jetGapFilter.cc
Go to the documentation of this file.
1 
8 
10 
12 
14 
15 
18 
22 
23 //
24 // constructors and destructor
25 //
27 {
28  inputTag_ = iConfig.getParameter<edm::InputTag> ("inputTag");
29  minEt_ = iConfig.getParameter<double> ("minEt");
30  minEta_ = iConfig.getParameter<double> ("minEta");
31 
32  m_theCaloJetToken = consumes<reco::CaloJetCollection>(inputTag_);
33 
34 }
35 
37 
38 void
42  desc.add<edm::InputTag>("inputTag",edm::InputTag("iterativeCone5CaloJets"));
43  desc.add<double>("minEt",90.0);
44  desc.add<double>("minEta",1.9);
45  descriptions.add("hlt2jetGapFilter",desc);
46 }
47 
48 // ------------ method called to produce the data ------------
49 bool
51 {
52  using namespace trigger;
53 
54  // The filter object
55  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
56 
58  iEvent.getByToken(m_theCaloJetToken,recocalojets);
59 
60  // look at all candidates, check cuts and add to filter object
61  int n(0);
62 
63  if(recocalojets->size() > 1){
64  // events with two or more jets
65 
66  double etjet1=0.;
67  double etjet2=0.;
68  double etajet1=0.;
69  double etajet2=0.;
70  int countjets =0;
71 
72  for (auto recocalojet = recocalojets->begin();
73  recocalojet<=(recocalojets->begin()+1); recocalojet++) {
74 
75  if(countjets==0) {
76  etjet1 = recocalojet->et();
77  etajet1 = recocalojet->eta();
78  }
79  if(countjets==1) {
80  etjet2 = recocalojet->et();
81  etajet2 = recocalojet->eta();
82  }
83  countjets++;
84  }
85 
86  if(etjet1>minEt_ && etjet2>minEt_ && (etajet1*etajet2)<0 && std::abs(etajet1)>minEta_ && std::abs(etajet2)>minEta_){
87  for (auto recocalojet = recocalojets->begin();
88  recocalojet<=(recocalojets->begin()+1); recocalojet++) {
89  reco::CaloJetRef ref(reco::CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet)));
90  filterproduct.addObject(TriggerJet,ref);
91  n++;
92  }
93  }
94 
95  } // events with two or more jets
96 
97 
98 
99  // filter decision
100  bool accept(n>=2);
101 
102  return accept;
103 }
edm::InputTag inputTag_
T getParameter(std::string const &) const
HLT2jetGapFilter(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
~HLT2jetGapFilter() override
edm::EDGetTokenT< reco::CaloJetCollection > m_theCaloJetToken
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45