CMS 3D CMS Logo

HLTAcoFilter.cc
Go to the documentation of this file.
1 
10 #include <string>
12 
14 
16 
18 
21 
25 
26 
27 //
28 // constructors and destructor
29 //
31 {
32  inputJetTag_ = iConfig.getParameter< edm::InputTag > ("inputJetTag");
33  inputMETTag_ = iConfig.getParameter< edm::InputTag > ("inputMETTag");
34  minDPhi_ = iConfig.getParameter<double> ("minDeltaPhi");
35  maxDPhi_ = iConfig.getParameter<double> ("maxDeltaPhi");
36  minEtjet1_ = iConfig.getParameter<double> ("minEtJet1");
37  minEtjet2_ = iConfig.getParameter<double> ("minEtJet2");
38  AcoString_ = iConfig.getParameter<std::string> ("Acoplanar");
39 
40  m_theJetToken = consumes<reco::CaloJetCollection>(inputJetTag_);
41  m_theMETToken = consumes<trigger::TriggerFilterObjectWithRefs>(inputMETTag_);
42 }
43 
45 
46 void
50  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("IterativeCone5CaloJets"));
51  desc.add<edm::InputTag>("inputMETTag",edm::InputTag("MET"));
52  desc.add<double>("minDeltaPhi",0.0);
53  desc.add<double>("maxDeltaPhi",2.0);
54  desc.add<double>("minEtJet1",20.0);
55  desc.add<double>("minEtJet2",20.0);
56  desc.add<std::string>("Acoplanar","Jet1Jet2");
57  descriptions.add("hltAcoFilter",desc);
58 }
59 
60 // ------------ method called to produce the data ------------
61 bool
63 {
64  using namespace std;
65  using namespace edm;
66  using namespace reco;
67  using namespace trigger;
68 
69  // The filter object
70  if (saveTags()) {
71  filterproduct.addCollectionTag(inputJetTag_);
72  filterproduct.addCollectionTag(inputMETTag_);
73  }
74 
75  Handle<CaloJetCollection> recocalojets;
76  iEvent.getByToken(m_theJetToken,recocalojets);
78  iEvent.getByToken(m_theMETToken,metcal);
79 
80  // look at all candidates, check cuts and add to filter object
81  int n(0);
82  int JetNum = recocalojets->size();
83 
84  // events with two or more jets
85 
86  double etjet1=0.;
87  double etjet2=0.;
88  double phijet1=0.;
89  double phijet2=0.;
90  //double etmiss=0.;
91  double phimiss=0.;
92 
93  VRcalomet vrefMET;
94  metcal->getObjects(TriggerMET,vrefMET);
95  CaloMETRef metRef=vrefMET.at(0);
96  //etmiss = vrefMET.at(0)->et();
97  phimiss = vrefMET.at(0)->phi();
98 
99  CaloJetRef ref1,ref2;
100 
101  if (JetNum>0) {
102  auto recocalojet = recocalojets->begin();
103 
104  etjet1 = recocalojet->et();
105  phijet1 = recocalojet->phi();
106  ref1 = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
107 
108  if(JetNum>1) {
109  recocalojet++;
110  etjet2 = recocalojet->et();
111  phijet2 = recocalojet->phi();
112  ref2 = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
113  }
114  double Dphi= -1.;
115  int JetSel = 0;
116 
117  if (AcoString_ == "Jet2Met") {
118  Dphi = std::abs(phimiss-phijet2);
119  if (JetNum>=2 && etjet1>minEtjet1_ && etjet2>minEtjet2_) {JetSel=1;}
120  }
121  if (AcoString_ == "Jet1Jet2") {
122  Dphi = std::abs(phijet1-phijet2);
123  if (JetNum>=2 && etjet1>minEtjet1_ && etjet2>minEtjet2_) {JetSel=1;}
124  }
125  if (AcoString_ == "Jet1Met") {
126  Dphi = std::abs(phimiss-phijet1);
127  if (JetNum>=1 && etjet1>minEtjet1_ ) {JetSel=1;}
128  }
129 
130 
131  if (Dphi>M_PI) {Dphi=2.0*M_PI-Dphi;}
132  if(JetSel>0 && Dphi>=minDPhi_ && Dphi<=maxDPhi_){
133 
134  if (AcoString_=="Jet2Met" || AcoString_=="Jet1Met") {filterproduct.addObject(TriggerMET,metRef);}
135  if (AcoString_=="Jet1Met" || AcoString_=="Jet1Jet2") {filterproduct.addObject(TriggerJet,ref1);}
136  if (AcoString_=="Jet2Met" || AcoString_=="Jet1Jet2") {filterproduct.addObject(TriggerJet,ref2);}
137  n++;
138  }
139 
140 
141  } // at least one jet
142 
143 
144  // filter decision
145  bool accept(n>=1);
146 
147  return accept;
148 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::CaloJetCollection > m_theJetToken
Definition: HLTAcoFilter.h:36
~HLTAcoFilter() override
std::string AcoString_
Definition: HLTAcoFilter.h:45
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
double minEtjet1_
Definition: HLTAcoFilter.h:41
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition: HLTAcoFilter.cc:62
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theMETToken
Definition: HLTAcoFilter.h:37
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
edm::InputTag inputJetTag_
Definition: HLTAcoFilter.h:39
edm::InputTag inputMETTag_
Definition: HLTAcoFilter.h:40
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
HLTAcoFilter(const edm::ParameterSet &)
Definition: HLTAcoFilter.cc:30
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double minEtjet2_
Definition: HLTAcoFilter.h:42
#define M_PI
double minDPhi_
Definition: HLTAcoFilter.h:43
edm::Ref< CaloJetCollection > CaloJetRef
edm references
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTAcoFilter.cc:47
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
double maxDPhi_
Definition: HLTAcoFilter.h:44