CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetSortedVBFFilter.cc
Go to the documentation of this file.
1 
12 #include <vector>
13 #include <string>
14 
24 
25 using namespace std;
26 //
27 // constructors and destructor//
28 //
29 template<typename T>
31 ,inputJets_ (iConfig.getParameter<edm::InputTag>("inputJets" ))
32 ,inputJetTags_(iConfig.getParameter<edm::InputTag>("inputJetTags"))
33 ,mqq_ (iConfig.getParameter<double> ("Mqq" ))
34 ,detaqq_ (iConfig.getParameter<double> ("Detaqq" ))
35 ,detabb_ (iConfig.getParameter<double> ("Detabb" ))
36 ,dphibb_ (iConfig.getParameter<double> ("Dphibb" ))
37 ,ptsqq_ (iConfig.getParameter<double> ("Ptsumqq" ))
38 ,ptsbb_ (iConfig.getParameter<double> ("Ptsumbb" ))
39 ,seta_ (iConfig.getParameter<double> ("Etaq1Etaq2" ))
40 ,njets_ (iConfig.getParameter<int> ("njets" ))
41 ,value_ (iConfig.getParameter<std::string> ("value" ))
42 ,triggerType_ (iConfig.getParameter<int> ("triggerType" ))
43 {
44  m_theJetsToken = consumes<std::vector<T>>(inputJets_);
45  m_theJetTagsToken = consumes<reco::JetTagCollection>(inputJetTags_);
46  if(njets_<4) {
47  edm::LogWarning("LowNJets")<< "njets="<<njets_<<" it must be >=4. Forced njets=4.";
48  njets_=4;
49  }
50 }
51 
52 
53  template<typename T>
55 { }
56 
57 template<typename T>
58 void
61  makeHLTFilterDescription(desc);
62  desc.add<edm::InputTag>("inputJets",edm::InputTag("hltJetCollection"));
63  desc.add<edm::InputTag>("inputJetTags",edm::InputTag(""));
64  desc.add<double>("Mqq",200);
65  desc.add<double>("Detaqq",2.5);
66  desc.add<double>("Detabb",10.);
67  desc.add<double>("Dphibb",10.);
68  desc.add<double>("Ptsumqq",0.);
69  desc.add<double>("Ptsumbb",0.);
70  desc.add<double>("Etaq1Etaq2",40.);
71  desc.add<std::string>("value","second");
72  desc.add<int>("triggerType",trigger::TriggerJet);
73  desc.add<int>("njets",4);
74  descriptions.add(defaultModuleLabel<HLTJetSortedVBFFilter<T>>(), desc);
75 }
76 
77 template<typename T> float HLTJetSortedVBFFilter<T>::findCSV(const typename std::vector<T>::const_iterator & jet, const reco::JetTagCollection & jetTags){
78  float minDr = 0.1;
79  float tmpCSV = -20 ;
80  for (reco::JetTagCollection::const_iterator jetb = jetTags.begin(); (jetb!=jetTags.end()); ++jetb) {
81  float tmpDr = reco::deltaR(*jet,*(jetb->first));
82 
83  if (tmpDr < minDr) {
84  minDr = tmpDr ;
85  tmpCSV= jetb->second;
86  }
87 
88  }
89  return tmpCSV;
90 
91 }
92 //
93 // member functions
94 //
95 
96 // ------------ method called to produce the data ------------
97 template<typename T>
98 bool
100 {
101  using namespace std;
102  using namespace edm;
103  using namespace reco;
104  using namespace trigger;
105 
106  typedef vector<T> TCollection;
107  typedef Ref<TCollection> TRef;
108 
109  bool accept(false);
110 
112  event.getByToken(m_theJetsToken,jets);
113  Handle<JetTagCollection> jetTags;
114 
115  if (saveTags()) filterproduct.addCollectionTag(inputJets_);
116  if (jets->size()<4) return false;
117 
118  const unsigned int nMax(njets_<jets->size()?njets_:jets->size());
119  vector<Jpair> sorted(nMax);
120  vector<TRef> jetRefs(nMax);
121  unsigned int nJetRefs=0;
122 
123  unsigned int nJet=0;
124  double value(0.0);
125 
127  if (inputJetTags_.encode()=="") {
128  for (typename TCollection::const_iterator jet=jets->begin(); (jet!=jets->end()&& nJet<nMax); ++jet) {
129  if (value_=="Pt") {
130  value=jet->pt();
131  } else if (value_=="Eta") {
132  value=jet->eta();
133  } else if (value_=="Phi") {
134  value=jet->phi();
135  } else {
136  value = 0.0;
137  }
138  sorted[nJet] = make_pair(value,nJet);
139  ++nJet;
140  }
141  sort(sorted.begin(),sorted.end(),comparator);
142  for (unsigned int i=0; i<nMax; ++i) {
143  jetRefs[i]=TRef(jets,sorted[i].second);
144  }
145  nJetRefs=nMax;
146  q1 = jetRefs[3]->p4();
147  b1 = jetRefs[2]->p4();
148  b2 = jetRefs[1]->p4();
149  q2 = jetRefs[0]->p4();
150  } else if(value_=="1BTagAndEta"){
151  event.getByToken(m_theJetTagsToken,jetTags);
152  vector<Jpair> sorted;
153  unsigned int b1_idx=-1;
154  float csv_max=-999;
155  for (typename TCollection::const_iterator jet=jets->begin(); (jet!=jets->end()&& nJet<nMax); ++jet) { //fill "sorted" and get the most b-tagged jet with higher CSV (b1)
156  value = findCSV(jet, *jetTags);
157  if(value>csv_max) {
158  csv_max=value;
159  b1_idx=nJet;
160  }
161  sorted.push_back(make_pair(jet->eta(),nJet));
162  nJet++;
163  // cout << "jetPt=" << jet->pt() << "\tjetEta=" << jet->eta() << "\tjetCSV=" << value << endl;
164  }
165  if(b1_idx>=sorted.size() || b1_idx<0) edm::LogError("OutOfRange")<< "b1 index out of range.";
166  sorted.erase(sorted.begin()+b1_idx); //remove the most b-tagged jet from "sorted"
167  sort(sorted.begin(),sorted.end(),comparator); //sort "sorted" by eta
168 
169  unsigned int q1_idx=sorted.front().second; //take the backward jet (q1)
170  unsigned int q2_idx=sorted.back().second; //take the forward jet (q2)
171 
172  unsigned int i=0;
173  while( (i==q1_idx) || (i==q2_idx) || (i==b1_idx) ) i++; //take jet with highest pT but q1,q2,b1 (q2)
174  unsigned int b2_idx=i;
175 
176  if(q1_idx<jets->size()) q1 = jets->at(q1_idx).p4(); else edm::LogWarning("Something wrong with q1");
177  if(q2_idx<jets->size()) q2 = jets->at(q2_idx).p4(); else edm::LogWarning("Something wrong with q2");
178  if(b1_idx<jets->size()) b1 = jets->at(b1_idx).p4(); else edm::LogWarning("Something wrong with b1");
179  if(b2_idx<jets->size()) b2 = jets->at(b2_idx).p4(); else edm::LogWarning("Something wrong with b2");
180 
181  jetRefs[0]= TRef(jets,b1_idx);
182  jetRefs[1]= TRef(jets,b2_idx);
183  jetRefs[2]= TRef(jets,q1_idx);
184  jetRefs[3]= TRef(jets,q2_idx);
185  nJetRefs=4;
186 
187  // cout<<"\tPathB: b1="<<b1.pt()<<" b2="<<b2.pt()<<" q1="<<q1.pt()<<" q2="<<q2.pt()<<endl;
188  } else if(value_=="2BTagAndPt"){
189  event.getByToken(m_theJetTagsToken,jetTags);
190  vector<Jpair> sorted;
191 
192  unsigned int b1_idx=-1;
193  unsigned int b2_idx=-1;
194  float csv1=-999;
195  float csv2=-999;
196  for (typename TCollection::const_iterator jet=jets->begin(); (jet!=jets->end()&& nJet<nMax); ++jet) { //fill "sorted" and get the two most b-tagged jets (b1,b2)
197  value = findCSV(jet, *jetTags);
198  if(value>csv1) {
199  csv2=csv1;
200  b2_idx=b1_idx;
201  csv1=value;
202  b1_idx=nJet;
203  }
204  else if(value>csv2){
205  csv2=value;
206  b2_idx=nJet;
207  }
208  sorted.push_back(make_pair(jet->eta(),nJet));
209  nJet++;
210  // cout << "jetPt=" << jet->pt() << "\tjetEta=" << jet->eta() << "\tjetCSV=" << value << endl;
211  }
212  sorted.erase(sorted.begin()+b1_idx); //remove b1 and b2 from sorted
213  sorted.erase(sorted.begin()+(b1_idx>b2_idx?b2_idx:b2_idx-1));
214 
215  unsigned int q1_idx=sorted.at(0).second; //get q1 and q2 as the jets with highest pT, but b1 and b2.
216  unsigned int q2_idx=sorted.at(1).second;
217 
218  if(q1_idx<jets->size()) q1 = jets->at(q1_idx).p4(); else edm::LogWarning("Something wrong with q1");
219  if(q2_idx<jets->size()) q2 = jets->at(q2_idx).p4(); else edm::LogWarning("Something wrong with q2");
220  if(b1_idx<jets->size()) b1 = jets->at(b1_idx).p4(); else edm::LogWarning("Something wrong with b1");
221  if(b2_idx<jets->size()) b2 = jets->at(b2_idx).p4(); else edm::LogWarning("Something wrong with b2");
222 
223  jetRefs[0]= TRef(jets,b1_idx);
224  jetRefs[1]= TRef(jets,b2_idx);
225  jetRefs[2]= TRef(jets,q1_idx);
226  jetRefs[3]= TRef(jets,q2_idx);
227  nJetRefs=4;
228 
229  // cout<<"\tPathA: b1="<<b1.pt()<<" b2="<<b2.pt()<<" q1="<<q1.pt()<<" q2="<<q2.pt()<<endl;
230  }
231  else {
232  event.getByToken(m_theJetTagsToken,jetTags);
233  for (typename TCollection::const_iterator jet=jets->begin(); (jet!=jets->end()&& nJet<nMax); ++jet) {
234 
235  if (value_=="second") {
236  value = findCSV(jet, *jetTags);
237  } else {
238  value = 0.0;
239  }
240  sorted[nJet] = make_pair(value,nJet);
241  ++nJet;
242  }
243  sort(sorted.begin(),sorted.end(),comparator);
244  for (unsigned int i=0; i<nMax; ++i) {
245  jetRefs[i]= TRef(jets,sorted[i].second);
246  }
247  nJetRefs=nMax;
248  b1 = jetRefs[3]->p4();
249  b2 = jetRefs[2]->p4();
250  q1 = jetRefs[1]->p4();
251  q2 = jetRefs[0]->p4();
252  }
253 
254  double mqq_bs = (q1+q2).M();
255  double deltaetaqq = std::abs(q1.Eta()-q2.Eta());
256  double deltaetabb = std::abs(b1.Eta()-b2.Eta());
257  double deltaphibb = std::abs(reco::deltaPhi(b1.Phi(),b2.Phi()));
258  double ptsqq_bs = (q1+q2).Pt();
259  double ptsbb_bs = (b1+b2).Pt();
260  double signeta = q1.Eta()*q2.Eta();
261 
262  if (
263  (mqq_bs > mqq_ ) &&
264  (deltaetaqq > detaqq_ ) &&
265  (deltaetabb < detabb_ ) &&
266  (deltaphibb < dphibb_ ) &&
267  (ptsqq_bs > ptsqq_ ) &&
268  (ptsbb_bs > ptsbb_ ) &&
269  (signeta < seta_ )
270  ) {
271  accept=true;
272  for (unsigned int i=0; i<nJetRefs; ++i) {
273  filterproduct.addObject(triggerType_,jetRefs[i]);
274  }
275  }
276 
277  return accept;
278 }
279 
edm::EDGetTokenT< reco::JetTagCollection > m_theJetTagsToken
std::string defaultModuleLabel()
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< std::vector< T > > m_theJetsToken
transient_vector_type::const_iterator const_iterator
const_iterator end() const
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
static float findCSV(const typename std::vector< T >::const_iterator &jet, const reco::JetTagCollection &jetTags)
double q2[4]
Definition: TauolaWrapper.h:88
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
U second(std::pair< T, U > const &p)
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
double q1[4]
Definition: TauolaWrapper.h:87
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLTJetSortedVBFFilter(const edm::ParameterSet &)
const_iterator begin() const
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
tuple size
Write out results.
math::PtEtaPhiELorentzVectorF LorentzVector
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override