CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLRJetCombObservables.cc
Go to the documentation of this file.
1 //
2 // Author: Jan Heyninck
3 // Created: Tue Apr 3 17:33:23 PDT 2007
4 //
5 //
12 
13 // constructor with path; default should not be used
15 : jetSourceToken_( jetSourceToken )
16 , genEvtToken_( iC.consumes<TtGenEvent>( edm::InputTag( "genEvt" ) ) )
17 {}
18 
19 
20 // destructor
22 
23 std::vector< TtSemiLRJetCombObservables::IntBoolPair >
25 {
26  bool debug=false;
27 
28  evtselectVarVal.clear();
29  evtselectVarMatch.clear();
30 
31  // Check whether the objects are matched:
32 // bool matchHadt = false; //commented out since gcc461 complained that this variable was set but unused
33 // bool matchLept = false; //commented out since gcc461 complained that this variable was set but unused
34 // bool matchHadW = false; //commented out since gcc461 complained that this variable was set but unused
35 // bool matchLepW = false; //commented out since gcc461 complained that this variable was set but unused
36  bool matchHadb = false;
37  bool matchLepb = false;
38  bool matchHadp = false;
39  bool matchHadq = false;
40  bool matchHadpq = false;
41  bool matchHadqp = false;
42  bool matchLepl = false;
43  bool matchLepn = false;
44 
45  if(debug) std::cout << "== start matching the objects " << std::endl;
46 
47  double drLepl=0, drLepn=0, drHadb=0, drLepb=0, drHadp=0, drHadq=0, drHadpq=0, drHadqp=0; // drHadt=0, drLept=0, drLepW=0, drHadW=0;
48 
50  iEvent.getByToken(genEvtToken_,genEvent);
51 
52  if (genEvent.failedToGet()) {
53  if(debug) std::cout << "== not found genEvent " << std::endl;
54  } else if (genEvent.isValid()) {
55  if(debug) std::cout << "== found genEvent " << std::endl;
56 
57  if (genEvent->isSemiLeptonic() && genEvent->numberOfBQuarks() == 2) {
58 
59  //if(debug) std::cout << "== genEvent->quarkFromAntiTop() " << genEvent->quarkFromAntiTop()->pt() << std::endl;
60  if(debug) std::cout << "== genEvent->isSemiLeptonic() && genEvent->numberOfBQuarks() == 2 " << std::endl;
61  if(debug) std::cout << "== solution.getDecay() == " <<solution.getDecay()<< std::endl;
62  if(debug) std::cout << "== solution.getRecLepm().pt() = " <<solution.getRecLepm().pt() << std::endl;
63  //if(debug) if(solution.getGenLepl() == 0) std::cout << "solution.getGenLepl() == NULL" << std::endl;
64  if(debug) std::cout << "== *(solution.getGenLept())" << solution.getGenLept()->pt() << std::endl;
65  if(debug) std::cout << "== *(solution.getGenLepl())" << solution.getGenLepl()->pt() << std::endl;
66  // std::cout << "Semilepton:\n";
67  // Match the lepton by deltaR
68  if (solution.getDecay() == "muon") drLepl = DeltaR<pat::Particle, reco::GenParticle>()(solution.getRecLepm(), *(solution.getGenLepl()));
69  if(debug) std::cout << "== matching lepton " << std::endl;
70  if (solution.getDecay() == "electron") drLepl = DeltaR<pat::Particle, reco::GenParticle>()(solution.getRecLepe(), *(solution.getGenLepl()));
71  matchLepl = (drLepl < 0.3);
72 
73  if(debug) std::cout << "== lepton is matched " << std::endl;
74  // Match the neutrino by deltaR
75  drLepn = DeltaR<pat::Particle, reco::GenParticle>()(solution.getRecLepn(), *(solution.getGenLepn()));
76  matchLepn = (drLepn < 0.3);
77 
78  // Match the hadronic b by deltaR
79  drHadb = DeltaR<pat::Jet, reco::GenParticle>()(solution.getRecHadb(), *(solution.getGenHadb()));
80  matchHadb = (drHadb < 0.3);
81 
82  // Match the hadronicleptonic b by deltaR
83  drLepb = DeltaR<pat::Jet, reco::GenParticle>()(solution.getRecLepb(), *(solution.getGenLepb()));
84  matchLepb = (drLepb < 0.3);
85 
86  // Match the hadronic p by deltaR
87  drHadp = DeltaR<pat::Jet, reco::GenParticle>()(solution.getRecHadp(), *(solution.getGenHadp()));
88  matchHadp = (drHadp < 0.3);
89 
90  // Match the hadronic pq by deltaR
91  drHadpq = DeltaR<pat::Jet, reco::GenParticle>()(solution.getRecHadp(), *(solution.getGenHadq()));
92  matchHadpq = (drHadpq < 0.3);
93 
94  // Match the hadronic q by deltaR
95  drHadq = DeltaR<pat::Jet, reco::GenParticle>()(solution.getRecHadq(), *(solution.getGenHadq()));
96  matchHadq = (drHadq < 0.3);
97 
98  // Match the hadronic qp by deltaR
99  drHadqp = DeltaR<pat::Jet, reco::GenParticle>()(solution.getRecHadq(), *(solution.getGenHadp()));
100  matchHadqp = (drHadqp < 0.3);
101 
102  //commented out since gcc461 complained that these variables were set but unused!
103  /*
104  // Match the hadronic W by deltaR
105  drHadW = DeltaR<reco::Particle, reco::GenParticle>()(solution.getRecHadW(), *(solution.getGenHadW()));
106  matchHadW = (drHadW < 0.3);
107 
108  // Match the leptonic W by deltaR
109  drLepW = DeltaR<reco::Particle, reco::GenParticle>()(solution.getRecLepW(), *(solution.getGenLepW()));
110  matchLepW = (drLepW < 0.3);
111 
112  // Match the hadronic t by deltaR
113  drHadt = DeltaR<reco::Particle, reco::GenParticle>()(solution.getRecHadt(), *(solution.getGenHadt()));
114  matchHadt = (drHadt < 0.3);
115 
116  // Match the leptonic t by deltaR
117  drLept = DeltaR<reco::Particle, reco::GenParticle>()(solution.getRecLept(), *(solution.getGenLept()));
118  matchLept = (drLept < 0.3);
119  */
120  }
121  } else {
122  if(debug) std::cout << "== nothing " << std::endl;
123  }
124 
125  if(debug) std::cout << "== objects matched" <<std::endl;
126 
128  iEvent.getByToken(jetSourceToken_, jets);
129 
130  if(debug) std::cout << "== start calculating observables" << std::endl;
131 
132 
133  //obs1 : pt(had top)
134  double AverageTop =((solution.getHadb().p4()+solution.getHadq().p4()+solution.getHadp().p4()).pt()+(solution.getLepb().p4()+solution.getHadq().p4()+solution.getHadp().p4()).pt()+(solution.getHadb().p4()+solution.getLepb().p4()+solution.getHadp().p4()).pt()+(solution.getHadb().p4()+solution.getHadq().p4()+solution.getLepb().p4()).pt())/4.;
135  double Obs1 = ((solution.getHadb().p4()+solution.getHadq().p4()+solution.getHadp().p4()).pt())/AverageTop;
136  evtselectVarVal.push_back(IntDblPair(1,Obs1));
137  evtselectVarMatch.push_back(IntBoolPair(1, ((matchHadq&&matchHadp)||(matchHadpq&&matchHadqp))&&matchHadb));
138 
139  if(debug) std::cout << "== observable 1 " << Obs1 << std::endl;
140 
141  //obs2 : (pt_b1 + pt_b2)/(sum jetpt)
142  double Obs2 = (solution.getHadb().pt()+solution.getLepb().pt())/(solution.getHadp().pt()+solution.getHadq().pt());
143  evtselectVarVal.push_back(IntDblPair(2,Obs2));
144  evtselectVarMatch.push_back(IntBoolPair(2,((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))&&matchHadb&&matchLepb));
145 
146  if(debug) std::cout << "== observable 2 " << Obs2 << std::endl;
147 
148  //obs3: delta R between lep b and lepton
149  double Obs3 = -999;
150  if (solution.getDecay() == "muon") Obs3 = ROOT::Math::VectorUtil::DeltaR( solution.getLepb().p4(),solution.getRecLepm().p4() );
151  if (solution.getDecay() == "electron") Obs3 = ROOT::Math::VectorUtil::DeltaR( solution.getLepb().p4(),solution.getRecLepe().p4() );
152  evtselectVarVal.push_back(IntDblPair(3,Obs3));
153  evtselectVarMatch.push_back(IntBoolPair(3,matchLepb&&matchLepl));
154 
155  if(debug) std::cout << "== observable 3 " << Obs3 << std::endl;
156 
157  //obs4 : del R ( had b, had W)
158  double Obs4 = ROOT::Math::VectorUtil::DeltaR( solution.getHadb().p4(), solution.getHadq().p4()+solution.getHadp().p4() );
159  evtselectVarVal.push_back(IntDblPair(4,Obs4));
160  evtselectVarMatch.push_back(IntBoolPair(4,matchHadb&&((matchHadp&&matchHadp)||(matchHadpq&&matchHadqp))));
161 
162  if(debug) std::cout << "== observable 4 " << Obs4 << std::endl;
163 
164  //obs5 : del R between light quarkssolution.getHadp().phi(
165  double Obs5 = ROOT::Math::VectorUtil::DeltaR( solution.getHadq().p4(),solution.getHadp().p4() );
166  evtselectVarVal.push_back(IntDblPair(5,Obs5));
167  evtselectVarMatch.push_back(IntBoolPair(5,(matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
168 
169  if(debug) std::cout << "== observable 5 " << Obs5 << std::endl;
170 
171  //obs6 : b-tagging information
172  double Obs6 = 0;
173  if ( fabs(solution.getHadb().bDiscriminator("trackCountingJetTags") +10) < 0.0001 || fabs(solution.getLepb().bDiscriminator("trackCountingJetTags") +10)< 0.0001 ){
174  Obs6 = -10.;
175  } else {
176  Obs6 = (solution.getHadb().bDiscriminator("trackCountingJetTags")+solution.getLepb().bDiscriminator("trackCountingJetTags"));
177  }
178  evtselectVarVal.push_back(IntDblPair(6,Obs6));
179  evtselectVarMatch.push_back(IntBoolPair(6,1));
180 
181  if(debug) std::cout << "== observable 6 " << Obs6 << std::endl;
182 
183  //obs7 : chi2 value of kinematical fit with W-mass constraint
184  double Obs7 =0;
185  if(solution.getProbChi2() <0){Obs7 = -0;} else { Obs7 = log10(solution.getProbChi2()+.00001);}
186  evtselectVarVal.push_back(IntDblPair(7,Obs7));
187  evtselectVarMatch.push_back(IntBoolPair(7,((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
188 
189  if(debug) std::cout << "== observable 7 " << Obs7 << std::endl;
190 
191  //obs8(=7+1)
192  double Obs8 = solution.getCalHadt().p4().pt();
193  evtselectVarVal.push_back(IntDblPair(8,Obs8));
194  evtselectVarMatch.push_back(IntBoolPair(8, matchHadb&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
195 
196  if(debug) std::cout << "== observable 8 " << Obs8 << std::endl;
197 
198  //obs9
199  double Obs9 = fabs(solution.getCalHadt().p4().eta());
200  evtselectVarVal.push_back(IntDblPair(9,Obs9));
201  evtselectVarMatch.push_back(IntBoolPair(9, matchHadb&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
202 
203  if(debug) std::cout << "== observable 9 " << Obs9 << std::endl;
204 
205  //obs10
206  double Obs10 = solution.getCalHadt().p4().theta();
207  evtselectVarVal.push_back(IntDblPair(10,Obs10));
208  evtselectVarMatch.push_back(IntBoolPair(10, matchHadb&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
209 
210  if(debug) std::cout << "== observable 10 " << Obs10 << std::endl;
211 
212  //obs11
213  double Obs11 = solution.getCalHadW().p4().pt();
214  evtselectVarVal.push_back(IntDblPair(11,Obs11));
215  evtselectVarMatch.push_back(IntBoolPair(11, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
216 
217  if(debug) std::cout << "== observable 11 " << Obs11 << std::endl;
218 
219  //obs12
220  double Obs12 = fabs(solution.getCalHadW().p4().eta());
221  evtselectVarVal.push_back(IntDblPair(12,Obs12));
222  evtselectVarMatch.push_back(IntBoolPair(12, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
223 
224  if(debug) std::cout << "== observable 12 " << Obs12 << std::endl;
225 
226  //obs13
227  double Obs13 = solution.getCalHadW().p4().theta();
228  evtselectVarVal.push_back(IntDblPair(13,Obs13));
229  evtselectVarMatch.push_back(IntBoolPair(13, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
230 
231  if(debug) std::cout << "== observable 13 " << Obs13 << std::endl;
232 
233  //obs14
234  double Obs14 = solution.getCalHadb().p4().pt();
235  evtselectVarVal.push_back(IntDblPair(14,Obs14));
236  evtselectVarMatch.push_back(IntBoolPair(14, matchHadb));
237 
238  if(debug) std::cout << "== observable 14 " << Obs14 << std::endl;
239 
240  //obs15
241  double Obs15 = fabs(solution.getCalHadb().p4().eta());
242  evtselectVarVal.push_back(IntDblPair(15,Obs15));
243  evtselectVarMatch.push_back(IntBoolPair(15, matchHadb));
244 
245  if(debug) std::cout << "== observable 15 " << Obs15 << std::endl;
246 
247  //obs16
248  double Obs16 = solution.getCalHadb().p4().theta();
249  evtselectVarVal.push_back(IntDblPair(16,Obs16));
250  evtselectVarMatch.push_back(IntBoolPair(16, matchHadb));
251 
252  if(debug) std::cout << "== observable 16 " << Obs16 << std::endl;
253 
254  //obs17
255  double Obs17 = solution.getCalHadp().p4().pt();
256  evtselectVarVal.push_back(IntDblPair(17,Obs17));
257  evtselectVarMatch.push_back(IntBoolPair(17, matchHadp));
258 
259  if(debug) std::cout << "== observable 17 " << Obs17 << std::endl;
260 
261  //obs18
262  double Obs18 = fabs(solution.getCalHadp().p4().eta());
263  evtselectVarVal.push_back(IntDblPair(18,Obs18));
264  evtselectVarMatch.push_back(IntBoolPair(18, matchHadp));
265 
266  if(debug) std::cout << "== observable 18 " << Obs18 << std::endl;
267 
268  //obs19
269  double Obs19 = solution.getCalHadp().p4().theta();
270  evtselectVarVal.push_back(IntDblPair(19,Obs19));
271  evtselectVarMatch.push_back(IntBoolPair(19, matchHadp));
272 
273  if(debug) std::cout << "== observable 19 " << Obs19 << std::endl;
274 
275  //obs20
276  double Obs20 = solution.getCalHadq().p4().pt();
277  evtselectVarVal.push_back(IntDblPair(20,Obs20));
278  evtselectVarMatch.push_back(IntBoolPair(20, matchHadq));
279 
280  if(debug) std::cout << "== observable 20 " << Obs20 << std::endl;
281 
282  //obs21
283  double Obs21 = fabs(solution.getCalHadq().p4().eta());
284  evtselectVarVal.push_back(IntDblPair(21,Obs21));
285  evtselectVarMatch.push_back(IntBoolPair(21, matchHadq));
286 
287  if(debug) std::cout << "== observable 21 " << Obs21 << std::endl;
288 
289  //obs22
290  double Obs22 = solution.getCalHadq().p4().theta();
291  evtselectVarVal.push_back(IntDblPair(22,Obs22));
292  evtselectVarMatch.push_back(IntBoolPair(22, matchHadq));
293 
294  if(debug) std::cout << "== observable 22 " << Obs22 << std::endl;
295 
296  //obs23
297  double Obs23 = solution.getCalLept().p4().pt();
298  evtselectVarVal.push_back(IntDblPair(23,Obs23));
299  evtselectVarMatch.push_back(IntBoolPair(23, matchLepl&&matchLepn&&matchLepb));
300 
301  if(debug) std::cout << "== observable 23 " << Obs23 << std::endl;
302 
303  //obs24
304  double Obs24 = fabs(solution.getCalLept().p4().eta());
305  evtselectVarVal.push_back(IntDblPair(24,Obs24));
306  evtselectVarMatch.push_back(IntBoolPair(24, matchLepl&&matchLepn&&matchLepb));
307 
308  if(debug) std::cout << "== observable 24 " << Obs24 << std::endl;
309 
310  //obs25
311  double Obs25 = solution.getCalLept().p4().theta();
312  evtselectVarVal.push_back(IntDblPair(25,Obs25));
313  evtselectVarMatch.push_back(IntBoolPair(25, matchLepl&&matchLepn&&matchLepb));
314 
315  if(debug) std::cout << "== observable 25 " << Obs25 << std::endl;
316 
317  //obs26
318  double Obs26 = solution.getRecLepW().p4().pt();
319  evtselectVarVal.push_back(IntDblPair(26,Obs26));
320  evtselectVarMatch.push_back(IntBoolPair(26, matchLepl&&matchLepn));
321 
322  if(debug) std::cout << "== observable 26 " << Obs26 << std::endl;
323 
324  //obs27
325  double Obs27 = fabs(solution.getRecLepW().p4().eta());
326  evtselectVarVal.push_back(IntDblPair(27,Obs27));
327  evtselectVarMatch.push_back(IntBoolPair(27, matchLepl&&matchLepn));
328 
329  if(debug) std::cout << "== observable 27 " << Obs27 << std::endl;
330 
331  //obs28
332  double Obs28 = solution.getRecLepW().p4().theta();
333  evtselectVarVal.push_back(IntDblPair(28,Obs28));
334  evtselectVarMatch.push_back(IntBoolPair(28, matchLepl&&matchLepn));
335 
336  if(debug) std::cout << "== observable 28 " << Obs28 << std::endl;
337 
338  //obs29
339  double Obs29 = solution.getCalLepb().p4().pt();
340  evtselectVarVal.push_back(IntDblPair(29,Obs29));
341  evtselectVarMatch.push_back(IntBoolPair(29, matchLepb));
342 
343  if(debug) std::cout << "== observable 29 " << Obs29 << std::endl;
344 
345  //obs30
346  double Obs30 = fabs(solution.getCalLepb().p4().eta());
347  evtselectVarVal.push_back(IntDblPair(30,Obs30));
348  evtselectVarMatch.push_back(IntBoolPair(30, matchLepb));
349 
350  if(debug) std::cout << "== observable 30 " << Obs30 << std::endl;
351 
352  //obs31
353  double Obs31 = solution.getCalLepb().p4().theta();
354  evtselectVarVal.push_back(IntDblPair(31,Obs31));
355  evtselectVarMatch.push_back(IntBoolPair(31, matchLepb));
356 
357  if(debug) std::cout << "== observable 31 " << Obs31 << std::endl;
358 
359  //obs32
360  double Obs32 = -999.;
361  if (solution.getDecay() == "muon") Obs32 = solution.getRecLepm().p4().pt();
362  if (solution.getDecay() == "electron") Obs32 = solution.getRecLepe().p4().pt();
363  evtselectVarVal.push_back(IntDblPair(32,Obs32));
364  evtselectVarMatch.push_back(IntBoolPair(32, matchLepl));
365 
366  if(debug) std::cout << "== observable 32 " << Obs32 << std::endl;
367 
368  //obs33
369  double Obs33 = -999.;
370  if (solution.getDecay() == "muon") Obs33 = fabs(solution.getRecLepm().p4().eta());
371  if (solution.getDecay() == "electron") Obs33 = fabs(solution.getRecLepe().p4().eta());
372  evtselectVarVal.push_back(IntDblPair(33,Obs33));
373  evtselectVarMatch.push_back(IntBoolPair(33, matchLepl));
374 
375  if(debug) std::cout << "== observable 33 " << Obs33 << std::endl;
376 
377  //obs34
378  double Obs34 = -999.;
379  if (solution.getDecay() == "muon") Obs34 = fabs(solution.getRecLepm().p4().theta());
380  if (solution.getDecay() == "electron") Obs34 = fabs(solution.getRecLepe().p4().theta());
381  evtselectVarVal.push_back(IntDblPair(34,Obs34));
382  evtselectVarMatch.push_back(IntBoolPair(34, matchLepl));
383 
384  if(debug) std::cout << "== observable 34 " << Obs34 << std::endl;
385 
386  //obs35
387  double Obs35 = solution.getFitLepn().p4().pt();
388  evtselectVarVal.push_back(IntDblPair(35,Obs35));
389  evtselectVarMatch.push_back(IntBoolPair(35, matchLepn));
390 
391  if(debug) std::cout << "== observable 35 " << Obs35 << std::endl;
392 
393  //obs36
394  double Obs36 = fabs(solution.getFitLepn().p4().eta());
395  evtselectVarVal.push_back(IntDblPair(36,Obs36));
396  evtselectVarMatch.push_back(IntBoolPair(36, matchLepn));
397 
398  if(debug) std::cout << "== observable 36 " << Obs36 << std::endl;
399 
400  //obs37
401  double Obs37 = solution.getFitLepn().p4().theta();
402  evtselectVarVal.push_back(IntDblPair(37,Obs37));
403  evtselectVarMatch.push_back(IntBoolPair(37, matchLepn));
404 
405  if(debug) std::cout << "== observable 37 " << Obs37 << std::endl;
406 
407  // 2 particle kinematics
408  //obs38
409  double Obs38 = fabs(solution.getCalHadW().p4().phi()- solution.getRecLepW().p4().phi());
410  if (Obs38 > 3.1415927) Obs38 = 2*3.1415927 - Obs31;
411  if (Obs38 < -3.1415927) Obs38 = -2*3.1415927 - Obs31;
412  evtselectVarVal.push_back(IntDblPair(38,Obs38));
413  evtselectVarMatch.push_back(IntBoolPair(38, matchLepl&&matchLepn&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
414 
415  if(debug) std::cout << "== observable 38 " << Obs38 << std::endl;
416 
417  //obs39
418  double Obs39 = fabs(solution.getCalHadW().p4().eta()- solution.getRecLepW().p4().eta());
419  evtselectVarVal.push_back(IntDblPair(39,Obs39));
420  evtselectVarMatch.push_back(IntBoolPair(39, matchLepl&&matchLepn&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
421 
422  if(debug) std::cout << "== observable 39 " << Obs39 << std::endl;
423 
424  //obs40
425  double Obs40 = fabs(solution.getCalHadW().p4().theta()- solution.getRecLepW().p4().theta());
426  evtselectVarVal.push_back(IntDblPair(40,Obs40));
427  evtselectVarMatch.push_back(IntBoolPair(40, matchLepl&&matchLepn&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
428 
429  if(debug) std::cout << "== observable 40 " << Obs40 << std::endl;
430 
431  //obs41
432  double Obs41 = ROOT::Math::VectorUtil::DeltaR(solution.getCalHadW().p4(), solution.getRecLepW().p4());
433  evtselectVarVal.push_back(IntDblPair(41,Obs41));
434  evtselectVarMatch.push_back(IntBoolPair(41, matchLepl&&matchLepn&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
435 
436  if(debug) std::cout << "== observable 41 " << Obs41 << std::endl;
437 
438  //obs42
439  double Obs42 = fabs(solution.getCalHadb().p4().phi()- solution.getCalLepb().p4().phi());
440  if (Obs42 > 3.1415927) Obs42 = 2*3.1415927 - Obs42;
441  if (Obs42 < -3.1415927) Obs42 = -2*3.1415927 - Obs42;
442  evtselectVarVal.push_back(IntDblPair(42,Obs42));
443  evtselectVarMatch.push_back(IntBoolPair(42, matchHadb&&matchLepb));
444 
445  if(debug) std::cout << "== observable 42 " << Obs42 << std::endl;
446 
447  //obs43
448  double Obs43 = fabs(solution.getCalHadb().p4().eta()- solution.getCalLepb().p4().eta());
449  evtselectVarVal.push_back(IntDblPair(43,Obs43));
450  evtselectVarMatch.push_back(IntBoolPair(43, matchHadb&&matchLepb));
451 
452  if(debug) std::cout << "== observable 43 " << Obs43 << std::endl;
453 
454  //obs44
455  double Obs44 = fabs(solution.getCalHadb().p4().theta()- solution.getCalLepb().p4().theta());
456  evtselectVarVal.push_back(IntDblPair(44,Obs44));
457  evtselectVarMatch.push_back(IntBoolPair(44, matchHadb&&matchLepb));
458 
459  if(debug) std::cout << "== observable 44 " << Obs44 << std::endl;
460 
461  //obs45
462  double Obs45 = ROOT::Math::VectorUtil::DeltaR(solution.getCalHadb().p4(), solution.getCalLepb().p4());
463  evtselectVarVal.push_back(IntDblPair(45,Obs45));
464  evtselectVarMatch.push_back(IntBoolPair(45, matchHadb&&matchLepb));
465 
466  if(debug) std::cout << "== observable 45 " << Obs45 << std::endl;
467 
468  //obs46
469  double Obs46 = fabs(solution.getCalHadb().p4().phi()- solution.getCalHadW().p4().phi());
470  if (Obs46 > 3.1415927) Obs46 = 2*3.1415927 - Obs46;
471  if (Obs46 < -3.1415927) Obs46 = -2*3.1415927 - Obs46;
472  evtselectVarVal.push_back(IntDblPair(46,Obs46));
473  evtselectVarMatch.push_back(IntBoolPair(46, matchHadb&&((matchHadq&&matchHadp)||(matchHadpq&&matchHadqp))));
474 
475  if(debug) std::cout << "== observable 46 " << Obs46 << std::endl;
476 
477  //obs47
478  double Obs47 = fabs(solution.getCalHadb().p4().eta()- solution.getCalHadW().p4().eta());
479  evtselectVarVal.push_back(IntDblPair(47,Obs47));
480  evtselectVarMatch.push_back(IntBoolPair(47, matchHadb&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
481 
482  if(debug) std::cout << "== observable 47 " << Obs47 << std::endl;
483 
484  //obs48
485  double Obs48 = fabs(solution.getCalHadb().p4().theta()- solution.getCalHadW().p4().theta());
486  evtselectVarVal.push_back(IntDblPair(48,Obs48));
487  evtselectVarMatch.push_back(IntBoolPair(48, matchHadb&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
488 
489  if(debug) std::cout << "== observable 48 " << Obs48 << std::endl;
490 
491  //obs49
492  double Obs49 = ROOT::Math::VectorUtil::DeltaR(solution.getCalHadb().p4(), solution.getCalHadW().p4());
493  evtselectVarVal.push_back(IntDblPair(49,Obs49));
494  evtselectVarMatch.push_back(IntBoolPair(49, matchHadb&&((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))));
495 
496  if(debug) std::cout << "== observable 49 " << Obs49 << std::endl;
497 
498  //obs50
499  double Obs50 = fabs(solution.getCalLepb().p4().phi()- solution.getRecLepW().p4().phi());
500  if (Obs50 > 3.1415927) Obs50 = 2*3.1415927 - Obs50;
501  if (Obs50 < -3.1415927) Obs50 = -2*3.1415927 - Obs50;
502  evtselectVarVal.push_back(IntDblPair(50,Obs50));
503  evtselectVarMatch.push_back(IntBoolPair(50, matchLepb&&matchLepn&&matchLepl));
504 
505  if(debug) std::cout << "== observable 50 " << Obs50 << std::endl;
506 
507  //obs51
508  double Obs51 = fabs(solution.getCalLepb().p4().eta()- solution.getRecLepW().p4().eta());
509  evtselectVarVal.push_back(IntDblPair(51,Obs51));
510  evtselectVarMatch.push_back(IntBoolPair(51, matchLepb&&matchLepn&&matchLepl));
511 
512  if(debug) std::cout << "== observable 51 " << Obs51 << std::endl;
513 
514  //obs52
515  double Obs52 = fabs(solution.getCalLepb().p4().theta()- solution.getRecLepW().p4().theta());
516  evtselectVarVal.push_back(IntDblPair(52,Obs52));
517  evtselectVarMatch.push_back(IntBoolPair(52, matchLepb&&matchLepn&&matchLepl));
518 
519  if(debug) std::cout << "== observable 52 " << Obs52 << std::endl;
520 
521  //obs53
522  double Obs53 = ROOT::Math::VectorUtil::DeltaR(solution.getCalLepb().p4(), solution.getRecLepW().p4());
523  evtselectVarVal.push_back(IntDblPair(53,Obs53));
524  evtselectVarMatch.push_back(IntBoolPair(53, matchLepb&&matchLepn&&matchLepl));
525 
526  if(debug) std::cout << "== observable 53 " << Obs53 << std::endl;
527 
528  //obs54
529  double Obs54 = fabs(solution.getCalHadp().p4().phi()- solution.getCalHadq().p4().phi());
530  if (Obs54 > 3.1415927) Obs54 = 2*3.1415927 - Obs54;
531  if (Obs54 < -3.1415927) Obs54 = -2*3.1415927 - Obs54;
532  evtselectVarVal.push_back(IntDblPair(54,Obs54));
533  evtselectVarMatch.push_back(IntBoolPair(54, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
534 
535  if(debug) std::cout << "== observable 54 " << Obs54 << std::endl;
536 
537  //obs55
538  double Obs55 = fabs(solution.getCalHadp().p4().eta()- solution.getCalHadq().p4().eta());
539  evtselectVarVal.push_back(IntDblPair(55,Obs55));
540  evtselectVarMatch.push_back(IntBoolPair(55, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
541 
542  if(debug) std::cout << "== observable 55 " << Obs55 << std::endl;
543 
544  //obs56
545  double Obs56 = fabs(solution.getCalHadp().p4().theta()- solution.getCalHadq().p4().theta());
546  evtselectVarVal.push_back(IntDblPair(56,Obs56));
547  evtselectVarMatch.push_back(IntBoolPair(56, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
548 
549  if(debug) std::cout << "== observable 56 " << Obs56 << std::endl;
550 
551  //obs57
552  double Obs57 = ROOT::Math::VectorUtil::DeltaR(solution.getCalHadp().p4(), solution.getCalHadq().p4());
553  evtselectVarVal.push_back(IntDblPair(57,Obs57));
554  evtselectVarMatch.push_back(IntBoolPair(57, (matchHadp&&matchHadq)||(matchHadpq&&matchHadqp)));
555 
556  if(debug) std::cout << "== observable 57 " << Obs57 << std::endl;
557 
558  //obs58
559  double Obs58 = -999.;
560  if (solution.getDecay() == "muon") Obs58 = fabs(solution.getRecLepm().p4().phi()- solution.getRecLepn().p4().phi());
561  if (solution.getDecay() == "electron") Obs58 = fabs(solution.getRecLepe().p4().phi()- solution.getRecLepn().p4().phi());
562  if (Obs58 > 3.1415927) Obs58 = 2*3.1415927 - Obs58;
563  if (Obs58 < -3.1415927) Obs58 = -2*3.1415927 - Obs58;
564  evtselectVarVal.push_back(IntDblPair(58,Obs58));
565  evtselectVarMatch.push_back(IntBoolPair(58, matchLepl&&matchLepn));
566 
567  if(debug) std::cout << "== observable 58 " << Obs58 << std::endl;
568 
569  //obs59
570  double Obs59 = -999.;
571  if (solution.getDecay() == "muon") Obs59 = fabs(solution.getRecLepm().p4().eta()- solution.getRecLepn().p4().eta());
572  if (solution.getDecay() == "electron") Obs59 = fabs(solution.getRecLepe().p4().eta()- solution.getRecLepn().p4().eta());
573  evtselectVarVal.push_back(IntDblPair(59,Obs59));
574  evtselectVarMatch.push_back(IntBoolPair(59, matchLepl&&matchLepn));
575 
576  if(debug) std::cout << "== observable 59 " << Obs59 << std::endl;
577 
578  //obs60
579  double Obs60 = -999.;
580  if (solution.getDecay() == "muon") Obs60 = fabs(solution.getRecLepm().p4().theta()- solution.getRecLepn().p4().theta());
581  if (solution.getDecay() == "electron") Obs60 = fabs(solution.getRecLepe().p4().theta()- solution.getRecLepn().p4().theta());
582  evtselectVarVal.push_back(IntDblPair(60,Obs60));
583  evtselectVarMatch.push_back(IntBoolPair(60, matchLepl&&matchLepn));
584 
585  if(debug) std::cout << "== observable 60 " << Obs60 << std::endl;
586 
587  //obs61
588  double Obs61 = -999.;
589  if (solution.getDecay() == "muon") Obs61 = ROOT::Math::VectorUtil::DeltaR(solution.getRecLepm().p4(), solution.getRecLepn().p4());
590  if (solution.getDecay() == "electron") Obs61 = ROOT::Math::VectorUtil::DeltaR(solution.getRecLepe().p4(), solution.getRecLepn().p4());
591  evtselectVarVal.push_back(IntDblPair(61,Obs61));
592  evtselectVarMatch.push_back(IntBoolPair(61, matchLepl&&matchLepn));
593 
594  if(debug) std::cout << "== observable 61 " << Obs61 << std::endl;
595 
596  // miscellaneous event variables
597 
598 
599  //obs62
600  double Obs62 = ((jets->size() > 4 && (*jets)[3].p4().Et() > 0.00001) ? (*jets)[4].p4().Et() / (*jets)[3].p4().Et() : 1.0);
601  //double Obs62 = 1;
602  evtselectVarVal.push_back(IntDblPair(62,Obs62));
603  evtselectVarMatch.push_back(IntBoolPair(62, ((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))&&matchHadb&&matchLepb));
604 
605  if(debug) std::cout << "== observable 62 " << Obs62 << std::endl;
606 
607  float calJetsSumEt = 0;
608  for (unsigned int i = 4; i < jets->size(); i++) {
609  calJetsSumEt += (*jets)[i].p4().Et();
610  }
611 
612  //obs63
613  double Obs63_den = (jets->size() > 4) ? ((*jets)[0].p4().Et()+(*jets)[1].p4().Et()+(*jets)[2].p4().Et()+(*jets)[3].p4().Et()+(*jets)[4].p4().Et()) : 0.0;
614  double Obs63 = (Obs63_den > 0.00001) ? calJetsSumEt / Obs63_den : 1.0;
615  //double Obs63 =1;
616  evtselectVarVal.push_back(IntDblPair(63,Obs63));
617  evtselectVarMatch.push_back(IntBoolPair(63, ((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))&&matchHadb&&matchLepb));
618 
619  if(debug) std::cout << "== observable 63 " << Obs63 << std::endl;
620 
621  //obs64
622  double Obs64 = solution.getProbChi2();
623  evtselectVarVal.push_back(IntDblPair(64,Obs64));
624  evtselectVarMatch.push_back(IntBoolPair(64, ((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))&&matchHadb&&matchLepb));
625 
626  if(debug) std::cout << "== observable 64 " << Obs64 << std::endl;
627 
628  //obs65
629  double Obs65 = solution.getFitHadt().p4().mass() - solution.getCalHadt().p4().mass();
630  evtselectVarVal.push_back(IntDblPair(65,Obs65));
631  evtselectVarMatch.push_back(IntBoolPair(65, ((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))&&matchHadb&&matchLepb));
632 
633  if(debug) std::cout << "== observable 65 " << Obs65 << std::endl;
634 
635  //obs66
636  double Obs66 = solution.getFitLept().p4().mass() - solution.getCalLept().p4().mass();
637  evtselectVarVal.push_back(IntDblPair(66,Obs66));
638  evtselectVarMatch.push_back(IntBoolPair(66, ((matchHadp&&matchHadq)||(matchHadpq&&matchHadqp))&&matchHadb&&matchLepb));
639 
640  if(debug) std::cout << "observables calculated" << std::endl;
641 
642  if (!matchOnly) solution.setLRJetCombObservables(evtselectVarVal);
643  return evtselectVarMatch;
644 }
int i
Definition: DBlmapReader.cc:9
reco::Particle getCalHadW() const
pat::Jet getLepb() const
pat::Jet getRecLepb() const
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:223
std::pair< unsigned int, double > IntDblPair
pat::Jet getHadq() const
TtSemiLRJetCombObservables(edm::ConsumesCollector &&iC, const edm::EDGetTokenT< std::vector< pat::Jet > > &jetSourceToken)
Definition: deltaR.h:79
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
pat::Jet getCalHadq() const
std::pair< unsigned int, bool > IntBoolPair
const reco::GenParticle * getGenLept() const
pat::Particle getFitLepn() const
std::vector< IntDblPair > evtselectVarVal
reco::Particle getFitHadt() const
std::string getDecay() const
float bDiscriminator(const std::string &theLabel) const
-— methods for accessing b-tagging info -—
const reco::GenParticle * getGenHadq() const
const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: Particle.h:72
const reco::GenParticle * getGenLepb() const
virtual double pt() const
transverse momentum
pat::Jet getHadb() const
edm::EDGetTokenT< TtGenEvent > genEvtToken_
pat::Electron getRecLepe() const
std::vector< IntBoolPair > operator()(TtSemiEvtSolution &, const edm::Event &iEvent, bool matchOnly=false)
const reco::GenParticle * getGenHadp() const
int iEvent
Definition: GenABIO.cc:230
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
pat::Jet getRecHadp() const
const reco::GenParticle * getGenHadb() const
vector< PseudoJet > jets
pat::Jet getHadp() const
const reco::GenParticle * getGenLepl() const
reco::Particle getCalLept() const
edm::EDGetTokenT< std::vector< pat::Jet > > jetSourceToken_
tuple genEvent
Definition: MCTruth.py:33
bool isValid() const
Definition: HandleBase.h:75
pat::Jet getRecHadb() const
void setLRJetCombObservables(const std::vector< std::pair< unsigned int, double > > &varval)
pat::Jet getCalHadb() const
bool failedToGet() const
Definition: HandleBase.h:79
#define debug
Definition: HDRShower.cc:19
pat::Jet getRecHadq() const
const reco::GenParticle * getGenLepn() const
pat::Muon getRecLepm() const
pat::Jet getCalLepb() const
tuple cout
Definition: gather_cfg.py:121
pat::MET getRecLepn() const
std::vector< IntBoolPair > evtselectVarMatch
reco::Particle getCalHadt() const
reco::Particle getRecLepW() const
double getProbChi2() const
genEvtToken_(mayConsume< TtGenEvent >(genEvt_))
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
reco::Particle getFitLept() const
pat::Jet getCalHadp() const