CMS 3D CMS Logo

HLTPhi2METFilter.cc
Go to the documentation of this file.
1 
9 
11 
13 
15 
18 
21 
22 //
23 // constructors and destructor
24 //
26 {
27  inputJetTag_ = iConfig.getParameter< edm::InputTag > ("inputJetTag");
28  inputMETTag_ = iConfig.getParameter< edm::InputTag > ("inputMETTag");
29  minDPhi_ = iConfig.getParameter<double> ("minDeltaPhi");
30  maxDPhi_ = iConfig.getParameter<double> ("maxDeltaPhi");
31  minEtjet1_= iConfig.getParameter<double> ("minEtJet1");
32  minEtjet2_= iConfig.getParameter<double> ("minEtJet2");
33  m_theJetToken = consumes<reco::CaloJetCollection>(inputJetTag_);
34  m_theMETToken = consumes<trigger::TriggerFilterObjectWithRefs>(inputMETTag_);
35 }
36 
38 
42  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("iterativeCone5CaloJets"));
43  desc.add<edm::InputTag>("inputMETTag",edm::InputTag("hlt1MET60"));
44  desc.add<double>("minDeltaPhi",0.377);
45  desc.add<double>("minEtJet2",60.);
46  descriptions.add("hltPhi2METFilter",desc);
47 }
48 
49 // ------------ method called to produce the data ------------
50 bool
52 {
53  using namespace std;
54  using namespace edm;
55  using namespace reco;
56  using namespace trigger;
57 
58  // The filter object
59  if (saveTags()) {
60  filterproduct.addCollectionTag(inputJetTag_);
61  filterproduct.addCollectionTag(inputMETTag_);
62  }
63 
64  Handle<CaloJetCollection> recocalojets;
65  iEvent.getByToken(m_theJetToken,recocalojets);
67  iEvent.getByToken(m_theMETToken,metcal);
68 
69  // look at all candidates, check cuts and add to filter object
70  int n(0);
71 
72  VRcalomet vrefMET;
73  metcal->getObjects(TriggerMET,vrefMET);
74  CaloMETRef metRef=vrefMET.at(0);
75  CaloJetRef ref1,ref2;
76 
77  if(recocalojets->size() > 1){
78  // events with two or more jets
79 
80  double etjet1=0.;
81  double etjet2=0.;
82  double phijet2=0.;
83  // double etmiss = vrefMET.at(0)->et();
84  double phimiss = vrefMET.at(0)->phi();
85  int countjets =0;
86 
87  for (auto recocalojet = recocalojets->begin();
88  recocalojet<=(recocalojets->begin()+1); recocalojet++) {
89 
90  if(countjets==0) {
91  etjet1 = recocalojet->et();
92  ref1 = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
93  }
94  if(countjets==1) {
95  etjet2 = recocalojet->et();
96  phijet2 = recocalojet->phi();
97  ref2 = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
98  }
99  countjets++;
100  }
101  double Dphi= std::abs(phimiss-phijet2);
102  if (Dphi>M_PI) Dphi=2.0*M_PI-Dphi;
103  if(etjet1>minEtjet1_ && etjet2>minEtjet2_ && Dphi>=minDPhi_ && Dphi<=maxDPhi_){
104  filterproduct.addObject(TriggerMET,metRef);
105  filterproduct.addObject(TriggerJet,ref1);
106  filterproduct.addObject(TriggerJet,ref2);
107  n++;
108  }
109 
110  } // events with two or more jets
111 
112 
113 
114  // filter decision
115  bool accept(n>=1);
116 
117  return accept;
118 }
edm::InputTag inputMETTag_
T getParameter(std::string const &) const
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HLTPhi2METFilter(const edm::ParameterSet &)
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>)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
int iEvent
Definition: GenABIO.cc:224
edm::InputTag inputJetTag_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< reco::CaloJetCollection > m_theJetToken
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define M_PI
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)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
std::vector< reco::CaloMETRef > VRcalomet
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theMETToken
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~HLTPhi2METFilter() override