CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Top_Decaykin.cc
Go to the documentation of this file.
1 //
2 //
3 // File: src/Top_Decaykin.cc
4 // Purpose: Calculate some kinematic quantities for ttbar events.
5 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
6 //
7 // CMSSW File : src/Top_Decaykin.cc
8 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
9 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
10 //
11 
12 
39 #include <cmath>
40 #include <algorithm>
41 #include <ostream>
42 
43 
44 using std::sqrt;
45 using std::abs;
46 using std::swap;
47 using std::ostream;
48 
49 
50 namespace hitfit {
51 
52 
53 namespace {
54 
55 
61 Fourvec leptons (const Lepjets_Event& ev)
62 //
63 // Purpose: Sum all leptons in EV.
64 //
65 // Inputs:
66 // ev - The event.
67 //
68 // Returns:
69 // The sum of all leptons in EV.
70 //
71 {
72  return (ev.sum (lepton_label) +
73  ev.sum (electron_label) +
74  ev.sum (muon_label));
75 }
76 
77 
78 } // unnamed namespace
79 
80 
82  double tmass,
83  double& nuz1, double& nuz2)
84 //
85 // Purpose: Solve for the neutrino longitudinal z-momentum that makes
86 // the leptonic top have mass TMASS.
87 //
88 // Inputs:
89 // ev - The event to solve.
90 // tmass - The desired top mass.
91 //
92 // Outputs:
93 // nuz1 - First solution (smaller absolute value).
94 // nuz2 - Second solution.
95 //
96 // Returns:
97 // True if there was a real solution. False if there were only
98 // imaginary solutions. (In that case, we just set the imaginary
99 // part to zero.)
100 //
101 {
102  bool discrim_flag = true;
103 
104  const Fourvec& vnu = ev.met ();
105  Fourvec cprime = leptons (ev) + ev.sum (lepb_label);
106  double alpha1 = tmass*tmass - cprime.m2();
107  double a = 2 * 4 * (cprime.z()*cprime.z() - cprime.e()*cprime.e());
108  double alpha = alpha1 + 2*(cprime.x()*vnu.x() + cprime.y()*vnu.y());
109  double b = 4 * alpha * cprime.z();
110  double c = alpha*alpha - 4 * cprime.e()*cprime.e() * vnu.vect().perp2();
111  double d = b*b - 2*a*c;
112  if (d < 0) {
113  discrim_flag = false;
114  d = 0;
115  }
116 
117  double dd = sqrt (d);
118  nuz1 = (-b + dd)/a;
119  nuz2 = (-b - dd)/a;
120  if (abs (nuz1) > abs (nuz2))
121  swap (nuz1, nuz2);
122 
123  return discrim_flag;
124 }
125 
126 
128  double tmass,
129  double& re_nuz1,
130  double& im_nuz1,
131  double& re_nuz2,
132  double& im_nuz2)
133 //
134 // Purpose: Solve for the neutrino longitudinal z-momentum that makes
135 // the leptonic top have mass TMASS, including the imaginary
136 // component of the z-component of the neutrino'somentum.
137 //
138 // Inputs:
139 // ev - The event to solve.
140 // tmass - The desired top mass.
141 //
142 // Outputs:
143 // re_nuz1 - Real component of the first solution.
144 // im_nuz1 - Imaginary component of the first solution (in the lower half of
145 // the complex plane).
146 // re_nuz2 - Real component of the second solution.
147 // im_nuz2 - Imaginary component of the second solution (in the upper half of
148 // the complex plane).
149 //
150 // Returns:
151 // True if there was a real solution. False if there were only
152 // complex solutions.
153 //
154 {
155  bool discrim_flag = true;
156 
157  const Fourvec& vnu = ev.met ();
158  Fourvec cprime = leptons (ev) + ev.sum (lepb_label);
159  double alpha1 = tmass*tmass - cprime.m2();
160  double a = 2 * 4 * (cprime.z()*cprime.z() - cprime.e()*cprime.e());
161  // Haryo's note: Here a is equivalent to '2a' in the quadratic
162  // equation ax^2 + bx + c = 0
163  double alpha = alpha1 + 2*(cprime.x()*vnu.x() + cprime.y()*vnu.y());
164  double b = 4 * alpha * cprime.z();
165  double c = alpha*alpha - 4 * cprime.e()*cprime.e() * vnu.vect().perp2();
166  double d = b*b - 2*a*c;
167  if (d < 0) {
168  discrim_flag = false;
169  }
170 
171  if (discrim_flag) {
172 
173  re_nuz1 = (-b + sqrt(d))/a ;
174  im_nuz1 = 0 ;
175  re_nuz2 = (-b - sqrt(d))/a ;
176  im_nuz2 = 0 ;
177  if (abs(re_nuz1) > abs(re_nuz2)) {
178  swap(re_nuz1,re_nuz2);
179  }
180 
181  } else {
182 
183  // Take absolute value of the imaginary component of nuz, in case
184  // a is negative, before multiplying by +1 or -1 to get the upper-half
185  // or lower-half imaginary value.
186  re_nuz1 = -b / a;
187  im_nuz1 = -fabs(sqrt(-d)/a );
188  re_nuz2 = -b / a;
189  im_nuz2 = fabs(sqrt(-d)/a );
190 
191 
192  }
193 
194 
195  return discrim_flag;
196 }
197 
198 
200  double wmass,
201  double& nuz1,
202  double& nuz2)
203 //
204 // Purpose: Solve for the neutrino longitudinal z-momentum that makes
205 // the leptonic W have mass WMASS.
206 //
207 // Inputs:
208 // ev - The event to solve.
209 // wmass - The desired W mass.
210 //
211 // Outputs:
212 // nuz1 - First solution (smaller absolute value).
213 // nuz2 - Second solution.
214 //
215 // Returns:
216 // True if there was a real solution. False if there were only
217 // imaginary solutions. (In that case, we just set the imaginary
218 // part to zero.)
219 //
220 {
221  bool discrim_flag = true;
222 
223  Fourvec vnu = ev.met();
224  Fourvec vlep = leptons (ev);
225 
226  double x = vlep.x()*vnu.x() + vlep.y()*vnu.y() + wmass*wmass/2;
227  double a = vlep.z()*vlep.z() - vlep.e()*vlep.e();
228  double b = 2*x*vlep.z();
229  double c = x*x - vnu.perp2() * vlep.e()*vlep.e();
230 
231  double d = b*b - 4*a*c;
232  if (d < 0) {
233  d = 0;
234  discrim_flag = false;
235  }
236 
237  nuz1 = (-b + sqrt (d))/2/a;
238  nuz2 = (-b - sqrt (d))/2/a;
239  if (abs (nuz1) > abs (nuz2))
240  swap (nuz1, nuz2);
241 
242  return discrim_flag;
243 }
244 
245 
247  double wmass,
248  double& re_nuz1,
249  double& im_nuz1,
250  double& re_nuz2,
251  double& im_nuz2)
252 //
253 // Purpose: Solve for the neutrino longitudinal z-momentum that makes
254 // the leptonic W have mass WMASS, including the imaginary
255 // component of the z-component of the neutrino'somentum.
256 //
257 // Inputs:
258 // ev - The event to solve.
259 // wmass - The desired W mass.
260 //
261 // Outputs:
262 // re_nuz1 - Real component of the first solution.
263 // im_nuz1 - Imaginary component of the first solution (in the lower half of
264 // the complex plane).
265 // re_nuz2 - Real component of the second solution.
266 // im_nuz2 - Imaginary component of the second solution (in the upper half of
267 // the complex plane).
268 //
269 // Returns:
270 // True if there was a real solution. False if there were only
271 // complex solutions.
272 //x
273 {
274  bool discrim_flag = true;
275 
276  Fourvec vnu = ev.met();
277  Fourvec vlep = leptons (ev);
278 
279  double x = vlep.x()*vnu.x() + vlep.y()*vnu.y() + wmass*wmass/2;
280  double a = vlep.z()*vlep.z() - vlep.e()*vlep.e();
281  double b = 2*x*vlep.z();
282  double c = x*x - vnu.perp2() * vlep.e()*vlep.e();
283 
284  double d = b*b - 4*a*c;
285  if (d < 0) {
286  discrim_flag = false;
287  }
288 
289  if (discrim_flag) {
290 
291  re_nuz1 = (-b + sqrt(d))/2/a ;
292  im_nuz1 = 0.0 ;
293  re_nuz2 = (-b - sqrt(d))/2/a ;
294  im_nuz2 = 0.0 ;
295  if (fabs(re_nuz1) > fabs(re_nuz2)) {
296  swap(re_nuz1,re_nuz2);
297  }
298 
299  } else {
300 
301  // Take absolute value of the imaginary component of nuz, in case
302  // a is negative, before multiplying by +1 or -1 to get the upper-half
303  // or lower-half imaginary value.
304 
305  re_nuz1 = -b /2/a ;
306  im_nuz1 = -fabs(sqrt(-d)/a);
307  re_nuz2 = -b /2/a ;
308  im_nuz2 = fabs(sqrt(-d)/a);
309 
310  }
311 
312  return discrim_flag;
313 }
314 
315 
317 //
318 // Purpose: Sum up the appropriate 4-vectors to find the hadronic W.
319 //
320 // Inputs:
321 // ev - The event.
322 //
323 // Returns:
324 // The hadronic W.
325 //
326 {
327  return (ev.sum (hadw1_label) + ev.sum (hadw2_label));
328 }
329 
330 
332 //
333 // Purpose:
334 //
335 // Inputs:
336 // ev - The event.
337 //
338 // Returns:
339 // The higher-pT hadronic jet from W
340 //
341 {
342  return ev.sum (hadw1_label);
343 }
344 
345 
347 //
348 // Purpose:
349 //
350 // Inputs:
351 // ev - The event.
352 //
353 // Returns:
354 // The lower-pT hadronic jet from W
355 //
356 //
357 {
358  return ev.sum (hadw2_label);
359 }
360 
361 
363 //
364 // Purpose: Sum up the appropriate 4-vectors to find the leptonic W.
365 //
366 // Inputs:
367 // ev - The event.
368 //
369 // Returns:
370 // The leptonic W.
371 //
372 {
373  return (leptons (ev) + ev.met ());
374 }
375 
376 
378 //
379 // Purpose: Sum up the appropriate 4-vectors to find the hadronic t.
380 //
381 // Inputs:
382 // ev - The event.
383 //
384 // Returns:
385 // The hadronic t.
386 //
387 {
388  return (ev.sum (hadb_label) + hadw (ev));
389 }
390 
391 
393 //
394 // Purpose: Sum up the appropriate 4-vectors to find the leptonic t.
395 //
396 // Inputs:
397 // ev - The event.
398 //
399 // Returns:
400 // The leptonic t.
401 //
402 {
403  return (ev.sum (lepb_label) + lepw (ev));
404 }
405 
406 
407 ostream& Top_Decaykin::dump_ev (std::ostream& s, const Lepjets_Event& ev)
408 //
409 // Purpose: Print kinematic information for EV.
410 //
411 // Inputs:
412 // s - The stream to which to write.
413 // ev - The event to dump.
414 //
415 // Returns:
416 // The stream S.
417 //
418 {
419  s << ev;
420  Fourvec p;
421 
422  p = lepw (ev);
423  s << "lepw " << p << " " << p.m() << "\n";
424  p = lept (ev);
425  s << "lept " << p << " " << p.m() << "\n";
426  p = hadw (ev);
427  s << "hadw " << p << " " << p.m() << "\n";
428  p = hadt (ev);
429  s << "hadt " << p << " " << p.m() << "\n";
430 
431  return s;
432 }
433 
434 
435 double Top_Decaykin::cos_theta_star(const Fourvec& fermion,
436  const Fourvec& W,
437  const Fourvec& top)
438 //
439 // Purpose: Calculate cos theta star in top decay
440 //
441 // Inputs:
442 // fermion - The four momentum of fermion from W
443 // W - The four momentum of W from top
444 // top - The four momentum of top
445 // Returns:
446 // cos theta star
447 //
448 {
449 
450  if (W.isLightlike() || W.isSpacelike()) {
451  return 100.0;
452  }
453 
454  CLHEP::HepBoost BoostWCM(W.findBoostToCM());
455 
456  CLHEP::Hep3Vector boost_v3fermion = BoostWCM(fermion).vect();
457  CLHEP::Hep3Vector boost_v3top = BoostWCM(top).vect();
458 
459  double costhetastar = boost_v3fermion.cosTheta(-boost_v3top);
460 
461  return costhetastar;
462 }
463 
464 
466 //
467 // Purpose: Calculate lepton cos theta star in top decay
468 //
469 // Inputs:
470 // ev - A lepton+jets event
471 // Returns:
472 // cos theta star of lepton
473 //
474 {
475 
476  return cos_theta_star(leptons(ev),
477  lepw(ev),
478  lept(ev));
479 
480 }
481 
482 
484 //
485 // Purpose: Calculate lepton cos theta star in top decay
486 //
487 // Inputs:
488 // ev - A lepton+jets event
489 // Returns:
490 // cos theta star of lepton
491 //
492 {
493 
494  return cos_theta_star(ev);
495 
496 }
497 
498 
500 //
501 // Purpose: Calculate absolute value of cos theta star of
502 // one of the hadronic W jet from hadronic top.
503 //
504 // Inputs:
505 // ev - A lepton+jets event
506 // Returns:
507 // absolute value of cos theta star
508 //
509 {
510 
511  return fabs(cos_theta_star(hadw1(ev),
512  hadw(ev),
513  hadt(ev)));
514 
515 }
516 
517 
518 } // namespace hitfit
519 
520 
static bool solve_nu_tmass(const Lepjets_Event &ev, double tmass, double &nuz1, double &nuz2)
Solve for the neutrino longitudinal momentum that makes the leptonic top have a certain value of mas...
Definition: Top_Decaykin.cc:81
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
float alpha
Definition: AMPTWrapper.h:95
Define three-vector and four-vector classes for the HitFit package, and supply a few additional opera...
static Fourvec hadw2(const Lepjets_Event &ev)
Return the hadronic boson jet which have lower .
bool ev
static Fourvec hadw1(const Lepjets_Event &ev)
Return the hadronic boson jet which have higher .
tuple d
Definition: ztail.py:151
static bool solve_nu(const Lepjets_Event &ev, double wmass, double &nuz1, double &nuz2)
Solve for the longitudinal momentum that makes the leptonic -boson to have a certain value of mass...
Represent a simple event consisting of lepton(s) and jet(s).
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:48
A class to hold functions to calculate kinematic quantities of interest in events.
static double cos_theta_star_lept(const Lepjets_Event &ev)
Calculate the lepton in top quark leptonic decay.
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:66
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Fourvec sum(int type) const
Return the sum of all objects&#39; four-momentum which have a particular type.
static double cos_theta_star_hadt(const Lepjets_Event &ev)
Calculate the hadronic in top quark leptonic decay. As there is no information on the weak isospin c...
Fourvec & met()
Return a reference to the missing transverse energy.
static Fourvec hadt(const Lepjets_Event &ev)
Sum up the appropriate four-momenta to find the hadronic top quark.
static double cos_theta_star(const Fourvec &fermion, const Fourvec &W, const Fourvec &top)
Calculate in top quark decay.
double b
Definition: hdecay.h:120
static Fourvec lept(const Lepjets_Event &ev)
Sum up the appropriate four-momenta to find the leptonic top quark.
double a
Definition: hdecay.h:121
static Fourvec lepw(const Lepjets_Event &ev)
Sum up the appropriate four-momenta to find the leptonic boson.
Definition: DDAxes.h:10
static Fourvec hadw(const Lepjets_Event &ev)
Sum up the appropriate four-momenta to find the hadronic boson.
static std::ostream & dump_ev(std::ostream &s, const Lepjets_Event &ev)
Print the kinematic information for an event.