CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTExclDiJetFilter.cc
Go to the documentation of this file.
1 
8 #include <cmath>
9 
21 
22 
23 //
24 // constructors and destructor
25 //
26 template<typename T>
28  HLTFilter(iConfig),
29  inputJetTag_ (iConfig.template getParameter<edm::InputTag> ("inputJetTag")),
30  caloTowerTag_(iConfig.template getParameter<edm::InputTag> ("caloTowerTag")),
31  minPtJet_ (iConfig.template getParameter<double> ("minPtJet")),
32  minHFe_ (iConfig.template getParameter<double> ("minHFe")),
33  HF_OR_ (iConfig.template getParameter<bool> ("HF_OR")),
34  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
35 {
36  m_theJetToken = consumes<std::vector<T>>(inputJetTag_);
37  m_theCaloTowerCollectionToken = consumes<CaloTowerCollection>(caloTowerTag_);
38  LogDebug("") << "HLTExclDiJetFilter: Input/minPtJet/minHFe/HF_OR/triggerType : "
39  << inputJetTag_.encode() << " "
40  << caloTowerTag_.encode() << " "
41  << minPtJet_ << " "
42  << minHFe_ << " "
43  << HF_OR_ << " "
44  << triggerType_;
45 }
46 
47 template<typename T>
49 
50 template<typename T>
51 void
54  makeHLTFilterDescription(desc);
55  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
56  desc.add<edm::InputTag>("caloTowerTag",edm::InputTag("hltTowerMakerForAll"));
57  desc.add<double>("minPtJet",30.0);
58  desc.add<double>("minHFe",50.0);
59  desc.add<bool>("HF_OR",false);
60  desc.add<int>("triggerType",trigger::TriggerJet);
61  descriptions.add(defaultModuleLabel<HLTExclDiJetFilter<T>>(), desc);
62 }
63 
64 // ------------ method called to produce the data ------------
65 template<typename T>
66 bool
68 {
69  using namespace std;
70  using namespace edm;
71  using namespace reco;
72  using namespace trigger;
73 
74  typedef vector<T> TCollection;
75  typedef Ref<TCollection> TRef;
76 
77  // The filter object
78  if (saveTags()) filterproduct.addCollectionTag(inputJetTag_);
79 
80  Handle<TCollection> recojets; //recojets can be any jet collections
81  iEvent.getByToken(m_theJetToken,recojets);
82 
83  // look at all candidates, check cuts and add to filter object
84  int n(0);
85 
86  double ptjet1=0., ptjet2=0.;
87  double phijet1=0., phijet2=0.;
88 
89  if(recojets->size() > 1){
90  // events with two or more jets
91 
92  int countjets =0;
93 
94  TRef JetRef1,JetRef2;
95 
96  typename TCollection::const_iterator recojet ( recojets->begin() );
97  for (;recojet<=(recojets->begin()+1); ++recojet) {
98  //
99  if(countjets==0) {
100  ptjet1 = recojet->pt();
101  phijet1 = recojet->phi();
102 
103  JetRef1 = TRef(recojets,distance(recojets->begin(),recojet));
104  }
105  //
106  if(countjets==1) {
107  ptjet2 = recojet->pt();
108  phijet2 = recojet->phi();
109 
110  JetRef2 = TRef(recojets,distance(recojets->begin(),recojet));
111  }
112  //
113  ++countjets;
114  }
115 
116  if(ptjet1>minPtJet_ && ptjet2>minPtJet_ ){
117  double Dphi=std::abs(phijet1-phijet2);
118  if(Dphi>M_PI) Dphi=2.0*M_PI-Dphi;
119  if(Dphi>0.5*M_PI) {
120  filterproduct.addObject(triggerType_,JetRef1);
121  filterproduct.addObject(triggerType_,JetRef2);
122  ++n;
123  }
124  } // pt(jet1,jet2) > minPtJet_
125 
126  } // events with two or more jets
127 
128  // calotowers
129  bool hf_accept=false;
130 
131  if(n>0) {
132  double ehfp(0.);
133  double ehfm(0.);
134 
136  iEvent.getByToken(m_theCaloTowerCollectionToken,o);
137 // if( o.isValid()) {
138  for( CaloTowerCollection::const_iterator cc = o->begin(); cc != o->end(); ++cc ) {
139  if(std::abs(cc->ieta())>28 && cc->energy()<4.0) continue;
140  if(cc->ieta()>28) ehfp+=cc->energy(); // HF+ energy
141  if(cc->ieta()<-28) ehfm+=cc->energy(); // HF- energy
142  }
143  // }
144 
145  bool hf_accept_and = (ehfp<minHFe_) && (ehfm<minHFe_);
146  bool hf_accept_or = (ehfp<minHFe_) || (ehfm<minHFe_);
147 
148  hf_accept = HF_OR_ ? hf_accept_or : hf_accept_and;
149 
150  } // n>0
151 
152 
154 
155 // filter decision
156  bool accept(n>0 && hf_accept);
157 
158  return accept;
159 }
#define LogDebug(id)
std::string defaultModuleLabel()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< CaloTower >::const_iterator const_iterator
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
std::string encode() const
Definition: InputTag.cc:164
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
edm::InputTag inputJetTag_
int iEvent
Definition: GenABIO.cc:230
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< std::vector< T > > m_theJetToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::InputTag caloTowerTag_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< CaloTowerCollection > m_theCaloTowerCollectionToken
#define M_PI
HLTExclDiJetFilter(const edm::ParameterSet &)
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) const override
def template
Definition: svgfig.py:520