CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTRapGapFilter.cc
Go to the documentation of this file.
1 
9 
11 
13 
15 
18 
21 
22 //
23 // constructors and destructor
24 //
26 {
27  inputTag_ = iConfig.getParameter< edm::InputTag > ("inputTag");
28  absEtaMin_ = iConfig.getParameter<double> ("minEta");
29  absEtaMax_ = iConfig.getParameter<double> ("maxEta");
30  caloThresh_ = iConfig.getParameter<double> ("caloThresh");
31  m_theJetToken = consumes<reco::CaloJetCollection>(inputTag_);
32 }
33 
35 
39  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("iterativeCone5CaloJets"));
40  desc.add<double>("minEta",3.0);
41  desc.add<double>("maxEta",5.0);
42  desc.add<double>("caloThresh",20.);
43  descriptions.add("hltRapGapFilter",desc);
44 }
45 
46 // ------------ method called to produce the data ------------
47 bool
49 {
50  using namespace reco;
51  using namespace trigger;
52 
53  // The filter object
54  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
55 
56  edm::Handle<CaloJetCollection> recocalojets;
57  iEvent.getByToken(m_theJetToken,recocalojets);
58 
59  // look at all candidates, check cuts and add to filter object
60  int n(0);
61 
62  //std::cout << "Found " << recocalojets->size() << " jets in this event" << std::endl;
63 
64  if(recocalojets->size() > 1){
65  // events with two or more jets
66 
67  double etjet=0.;
68  double etajet=0.;
69  double sumets=0.;
70  int countjets =0;
71 
72  for (CaloJetCollection::const_iterator recocalojet = recocalojets->begin();
73  recocalojet!=(recocalojets->end()); recocalojet++) {
74 
75  etjet = recocalojet->energy();
76  etajet = recocalojet->eta();
77 
78  if(std::abs(etajet) > absEtaMin_ && std::abs(etajet) < absEtaMax_)
79  {
80  sumets += etjet;
81  //std::cout << "Adding jet with eta = " << etajet << ", and e = "
82  // << etjet << std::endl;
83  }
84  countjets++;
85  }
86 
87  //std::cout << "Sum jet energy = " << sumets << std::endl;
88  if(sumets<=caloThresh_){
89  //std::cout << "Passed filter!" << std::endl;
90  for (CaloJetCollection::const_iterator recocalojet = recocalojets->begin();
91  recocalojet!=(recocalojets->end()); recocalojet++) {
92  CaloJetRef ref(CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet)));
93  filterproduct.addObject(TriggerJet,ref);
94  n++;
95  }
96  }
97 
98  } // events with two or more jets
99 
100  // filter decision
101  bool accept(n>0);
102 
103  return accept;
104 }
T getParameter(std::string const &) const
edm::InputTag inputTag_
HLTRapGapFilter(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
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:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::Ref< CaloJetCollection > CaloJetRef
edm references
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::CaloJetCollection > m_theJetToken