CMS 3D CMS Logo

Static Public Member Functions

hitfit::Top_Decaykin Class Reference

A class to hold functions to calculate kinematic quantities of interest in $ t\bar{t} \to \ell + 4 \mathrm{jets} $ events. This class has no state, only static member functions. More...

#include <Top_Decaykin.h>

List of all members.

Static Public Member Functions

static double cos_theta_star (const Fourvec fermion, const Fourvec W, const Fourvec top)
 Calculate $ \cos \theta^{*} $ in top quark decay.
static double cos_theta_star (const Lepjets_Event &ev)
 Calculate the lepton $ \cos \theta^{*} $ in top quark leptonic decay.
static double cos_theta_star_hadt (const Lepjets_Event &ev)
 Calculate the hadronic $ \cos \theta^{*} $ in top quark leptonic decay. As there is no information on the weak isospin component of the fermion, the absolute value of $ \cos \theta^{*} $ will be returned (the solutions for up-type and down-type fermions will differ only in sign but not in magnitude).
static double cos_theta_star_lept (const Lepjets_Event &ev)
 Calculate the lepton $ \cos \theta^{*} $ in top quark leptonic decay.
static std::ostream & dump_ev (std::ostream &s, const Lepjets_Event &ev)
 Print the kinematic information for an event.
static Fourvec hadt (const Lepjets_Event &ev)
 Sum up the appropriate four-momenta to find the hadronic top quark.
static Fourvec hadw (const Lepjets_Event &ev)
 Sum up the appropriate four-momenta to find the hadronic $ W- $ boson.
static Fourvec hadw1 (const Lepjets_Event &ev)
 Return the hadronic $ W- $ boson jet which have higher $ p_{T} $ .
static Fourvec hadw2 (const Lepjets_Event &ev)
 Return the hadronic $ W- $ boson jet which have lower $ p_{T} $ .
static Fourvec lept (const Lepjets_Event &ev)
 Sum up the appropriate four-momenta to find the leptonic top quark.
static Fourvec lepw (const Lepjets_Event &ev)
 Sum up the appropriate four-momenta to find the leptonic $ W- $ boson.
static bool solve_nu (const Lepjets_Event &ev, double wmass, double &re_nuz1, double &im_nuz1, double &re_nuz2, double &im_nuz2)
 Solve for the longitudinal $ z- $ momentum that makes the leptonic $ W $ -boson to have a certain value of mass. The complex component of the solutions are also given. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of real solutions, the first solution is the one which have smaller absolute value. In case of imaginary solutions (which are complex conjugate of each other), the first solution is the one which have imaginary component in the lower half of the complex plane, i.e., the one which have negative imaginary component).
static bool solve_nu (const Lepjets_Event &ev, double wmass, double &nuz1, double &nuz2)
 Solve for the longitudinal $ z- $ momentum that makes the leptonic $ W $ -boson to have a certain value of mass. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of complex solutions, the real components of the solutions are given.
static bool solve_nu_tmass (const Lepjets_Event &ev, double tmass, double &re_nuz1, double &im_nuz1, double &re_nuz2, double &im_nuz2)
 Solve for the neutrino longitudinal $ z- $ momentum that makes the leptonic top have a certain value of mass. The complex component of the solutions are also given. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of real solutions, the first solution is the one which have smaller absolute value. In case of imaginary solutions (which are complex conjugate of each other), the first solution is the one which have imaginary component in the lower half of the complex plane, i.e., the one which have negative imaginary component).
static bool solve_nu_tmass (const Lepjets_Event &ev, double tmass, double &nuz1, double &nuz2)
 Solve for the neutrino longitudinal $ z- $ momentum that makes the leptonic top have a certain value of mass. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of complex solutions, the real components of the solutions are given.

Detailed Description

A class to hold functions to calculate kinematic quantities of interest in $ t\bar{t} \to \ell + 4 \mathrm{jets} $ events. This class has no state, only static member functions.

Definition at line 56 of file Top_Decaykin.h.


Member Function Documentation

double hitfit::Top_Decaykin::cos_theta_star ( const Fourvec  fermion,
const Fourvec  W,
const Fourvec  top 
) [static]

Calculate $ \cos \theta^{*} $ in top quark decay.

Parameters:
fermionThe four-momentum of fermion from $ W- $ boson from top decay.
WThe four-momentum of $ W $ boson from top decay.
topThe four-momentum of top.

Definition at line 436 of file Top_Decaykin.cc.

Referenced by cos_theta_star(), cos_theta_star_hadt(), and cos_theta_star_lept().

{

    if (W.isLightlike() || W.isSpacelike()) {
        return 100.0;
    }

    CLHEP::HepBoost BoostWCM(W.findBoostToCM());

    CLHEP::Hep3Vector boost_v3fermion       = BoostWCM(fermion).vect();
    CLHEP::Hep3Vector boost_v3top           = BoostWCM(top).vect();

    double costhetastar = boost_v3fermion.cosTheta(-boost_v3top);

    return costhetastar;
}
double hitfit::Top_Decaykin::cos_theta_star ( const Lepjets_Event ev) [static]

Calculate the lepton $ \cos \theta^{*} $ in top quark leptonic decay.

Parameters:
evThe event to solve.

Definition at line 466 of file Top_Decaykin.cc.

References cos_theta_star(), lept(), EgammaValidation_Wenu_cff::leptons, and lepw().

{

    return cos_theta_star(leptons(ev),
                          lepw(ev),
                          lept(ev));

}
double hitfit::Top_Decaykin::cos_theta_star_hadt ( const Lepjets_Event ev) [static]

Calculate the hadronic $ \cos \theta^{*} $ in top quark leptonic decay. As there is no information on the weak isospin component of the fermion, the absolute value of $ \cos \theta^{*} $ will be returned (the solutions for up-type and down-type fermions will differ only in sign but not in magnitude).

Parameters:
evThe event to solve.

Definition at line 500 of file Top_Decaykin.cc.

References cos_theta_star(), hadt(), hadw(), and hadw1().

{

    return fabs(cos_theta_star(hadw1(ev),
                               hadw(ev),
                               hadt(ev)));

}
double hitfit::Top_Decaykin::cos_theta_star_lept ( const Lepjets_Event ev) [static]

Calculate the lepton $ \cos \theta^{*} $ in top quark leptonic decay.

Parameters:
evThe event to solve.

Definition at line 484 of file Top_Decaykin.cc.

References cos_theta_star().

{

    return cos_theta_star(ev);

}
ostream & hitfit::Top_Decaykin::dump_ev ( std::ostream &  s,
const Lepjets_Event ev 
) [static]

Print the kinematic information for an event.

Parameters:
sThe stream of which to write.
evThe event to be printed.

Definition at line 408 of file Top_Decaykin.cc.

References hadt(), hadw(), lept(), lepw(), AlCaHLTBitMon_ParallelJobs::p, and alignCSCRings::s.

Referenced by hitfit::Top_Fit::fit_one_perm().

{
  s << ev;
  Fourvec p;

  p = lepw (ev);
  s << "lepw " << p << " " << p.m() << "\n";
  p = lept (ev);
  s << "lept " << p << " " << p.m() << "\n";
  p = hadw (ev);
  s << "hadw " << p << " " << p.m() << "\n";
  p = hadt (ev);
  s << "hadt " << p << " " << p.m() << "\n";

  return s;
}
Fourvec hitfit::Top_Decaykin::hadt ( const Lepjets_Event ev) [static]

Sum up the appropriate four-momenta to find the hadronic top quark.

Parameters:
evThe event.

Definition at line 378 of file Top_Decaykin.cc.

References hitfit::hadb_label, hadw(), and hitfit::Lepjets_Event::sum().

Referenced by cos_theta_star_hadt(), dump_ev(), and hitfit::Top_Fit::fit_one_perm().

{
  return (ev.sum (hadb_label) + hadw (ev));
}
Fourvec hitfit::Top_Decaykin::hadw ( const Lepjets_Event ev) [static]

Sum up the appropriate four-momenta to find the hadronic $ W- $ boson.

Parameters:
evThe event.

Definition at line 317 of file Top_Decaykin.cc.

References hitfit::hadw1_label, hitfit::hadw2_label, and hitfit::Lepjets_Event::sum().

Referenced by cos_theta_star_hadt(), dump_ev(), hitfit::Top_Fit::fit_one_perm(), and hadt().

{
  return (ev.sum (hadw1_label) + ev.sum (hadw2_label));
}
Fourvec hitfit::Top_Decaykin::hadw1 ( const Lepjets_Event ev) [static]

Return the hadronic $ W- $ boson jet which have higher $ p_{T} $ .

Parameters:
evThe event.

Definition at line 332 of file Top_Decaykin.cc.

References hitfit::hadw1_label, and hitfit::Lepjets_Event::sum().

Referenced by cos_theta_star_hadt().

{
  return ev.sum (hadw1_label);
}
Fourvec hitfit::Top_Decaykin::hadw2 ( const Lepjets_Event ev) [static]

Return the hadronic $ W- $ boson jet which have lower $ p_{T} $ .

Parameters:
evThe event.

Definition at line 347 of file Top_Decaykin.cc.

References hitfit::hadw2_label, and hitfit::Lepjets_Event::sum().

{
  return ev.sum (hadw2_label);
}
Fourvec hitfit::Top_Decaykin::lept ( const Lepjets_Event ev) [static]

Sum up the appropriate four-momenta to find the leptonic top quark.

Parameters:
evThe event.

Definition at line 393 of file Top_Decaykin.cc.

References hitfit::lepb_label, lepw(), and hitfit::Lepjets_Event::sum().

Referenced by cos_theta_star(), dump_ev(), and hitfit::Top_Fit::fit_one_perm().

{
  return (ev.sum (lepb_label) + lepw (ev));
}
Fourvec hitfit::Top_Decaykin::lepw ( const Lepjets_Event ev) [static]

Sum up the appropriate four-momenta to find the leptonic $ W- $ boson.

Parameters:
evThe event.

Definition at line 363 of file Top_Decaykin.cc.

References EgammaValidation_Wenu_cff::leptons, and hitfit::Lepjets_Event::met().

Referenced by cos_theta_star(), dump_ev(), and lept().

{
  return (leptons (ev) + ev.met ());
}
bool hitfit::Top_Decaykin::solve_nu ( const Lepjets_Event ev,
double  wmass,
double &  re_nuz1,
double &  im_nuz1,
double &  re_nuz2,
double &  im_nuz2 
) [static]

Solve for the longitudinal $ z- $ momentum that makes the leptonic $ W $ -boson to have a certain value of mass. The complex component of the solutions are also given. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of real solutions, the first solution is the one which have smaller absolute value. In case of imaginary solutions (which are complex conjugate of each other), the first solution is the one which have imaginary component in the lower half of the complex plane, i.e., the one which have negative imaginary component).

Parameters:
evInput: The event to solve.
wmassInput: The desired mass of the $ W- $ boson in GeV.
re_nuz1Output: Real component of the first solution.
im_nuz1Output: Imaginary component of the first solution.
re_nuz2Output: Real component of the second solution.
im_nuz2Output: Imaginary component of the second solution.

Definition at line 247 of file Top_Decaykin.cc.

References a, b, trackerHits::c, EgammaValidation_Wenu_cff::leptons, hitfit::Lepjets_Event::met(), mathSSE::sqrt(), swap(), and x.

{
  bool discrim_flag = true;

  Fourvec vnu  = ev.met();
  Fourvec vlep = leptons (ev);

  double x = vlep.x()*vnu.x() + vlep.y()*vnu.y() + wmass*wmass/2;
  double a = vlep.z()*vlep.z() - vlep.e()*vlep.e();
  double b = 2*x*vlep.z();
  double c = x*x - vnu.perp2() * vlep.e()*vlep.e();

  double d = b*b - 4*a*c;
  if (d < 0) {
    discrim_flag = false;
  }

  if (discrim_flag) {

      re_nuz1 = (-b + sqrt(d))/2/a ;
      im_nuz1 = 0.0 ;
      re_nuz2 = (-b - sqrt(d))/2/a ;
      im_nuz2 = 0.0 ;
      if (fabs(re_nuz1) > fabs(re_nuz2)) {
          swap(re_nuz1,re_nuz2);
      }

  } else {

      // Take absolute value of the imaginary component of nuz, in case
      // a is negative, before multiplying by +1 or -1 to get the upper-half
      // or lower-half imaginary value.

      re_nuz1 = -b /2/a ;
      im_nuz1 = -fabs(sqrt(-d)/a);
      re_nuz2 = -b /2/a ;
      im_nuz2 =  fabs(sqrt(-d)/a);

  }

  return discrim_flag;
}
bool hitfit::Top_Decaykin::solve_nu ( const Lepjets_Event ev,
double  wmass,
double &  nuz1,
double &  nuz2 
) [static]

Solve for the longitudinal $ z- $ momentum that makes the leptonic $ W $ -boson to have a certain value of mass. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of complex solutions, the real components of the solutions are given.

Parameters:
evInput: The event to solve.
wmassInput: The desired mass of the $ W- $ boson in GeV.
nuz1Output: First solution (smaller absolute value).
nuz2Output: Second solution.

Definition at line 200 of file Top_Decaykin.cc.

References a, abs, b, trackerHits::c, EgammaValidation_Wenu_cff::leptons, hitfit::Lepjets_Event::met(), mathSSE::sqrt(), swap(), and x.

Referenced by hitfit::Top_Fit::fit_one_perm().

{
  bool discrim_flag = true;

  Fourvec vnu  = ev.met();
  Fourvec vlep = leptons (ev);

  double x = vlep.x()*vnu.x() + vlep.y()*vnu.y() + wmass*wmass/2;
  double a = vlep.z()*vlep.z() - vlep.e()*vlep.e();
  double b = 2*x*vlep.z();
  double c = x*x - vnu.perp2() * vlep.e()*vlep.e();

  double d = b*b - 4*a*c;
  if (d < 0) {
    d = 0;
    discrim_flag = false;
  }

  nuz1 = (-b + sqrt (d))/2/a;
  nuz2 = (-b - sqrt (d))/2/a;
  if (abs (nuz1) > abs (nuz2))
    swap (nuz1, nuz2);

  return discrim_flag;
}
bool hitfit::Top_Decaykin::solve_nu_tmass ( const Lepjets_Event ev,
double  tmass,
double &  re_nuz1,
double &  im_nuz1,
double &  re_nuz2,
double &  im_nuz2 
) [static]

Solve for the neutrino longitudinal $ z- $ momentum that makes the leptonic top have a certain value of mass. The complex component of the solutions are also given. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of real solutions, the first solution is the one which have smaller absolute value. In case of imaginary solutions (which are complex conjugate of each other), the first solution is the one which have imaginary component in the lower half of the complex plane, i.e., the one which have negative imaginary component).

Parameters:
evInput:The event to solve.
tmassInput: The desired value of top quark mass in GeV.
re_nuz1Output: Real component of the first solution.
im_nuz1Output: Imaginary component of the first solution.
re_nuz2Output: Real component of the second solution.
im_nuz2Output: Imaginary component of the second solution.

Definition at line 128 of file Top_Decaykin.cc.

References a, abs, alpha, b, trackerHits::c, hitfit::lepb_label, EgammaValidation_Wenu_cff::leptons, hitfit::Lepjets_Event::met(), mathSSE::sqrt(), hitfit::Lepjets_Event::sum(), and swap().

{
  bool discrim_flag = true;

  const Fourvec& vnu = ev.met ();
  Fourvec cprime = leptons (ev) + ev.sum (lepb_label);
  double alpha1 = tmass*tmass - cprime.m2();
  double a = 2 * 4 * (cprime.z()*cprime.z() - cprime.e()*cprime.e());
  // Haryo's note: Here a is equivalent to '2a' in the quadratic
  // equation ax^2 + bx + c = 0
  double alpha = alpha1 + 2*(cprime.x()*vnu.x() + cprime.y()*vnu.y());
  double b = 4 * alpha * cprime.z();
  double c = alpha*alpha - 4 * cprime.e()*cprime.e() * vnu.vect().perp2();
  double d = b*b - 2*a*c;
  if (d < 0) {
    discrim_flag = false;
  }

  if (discrim_flag) {

      re_nuz1 = (-b + sqrt(d))/a ;
      im_nuz1 = 0 ;
      re_nuz2 = (-b - sqrt(d))/a ;
      im_nuz2 = 0 ;
      if (abs(re_nuz1) > abs(re_nuz2)) {
          swap(re_nuz1,re_nuz2);
      }

  } else {

      // Take absolute value of the imaginary component of nuz, in case
      // a is negative, before multiplying by +1 or -1 to get the upper-half
      // or lower-half imaginary value.
      re_nuz1   = -b / a;
      im_nuz1   = -fabs(sqrt(-d)/a );
      re_nuz2   = -b / a;
      im_nuz2   =  fabs(sqrt(-d)/a );


  }


  return discrim_flag;
}
bool hitfit::Top_Decaykin::solve_nu_tmass ( const Lepjets_Event ev,
double  tmass,
double &  nuz1,
double &  nuz2 
) [static]

Solve for the neutrino longitudinal $ z- $ momentum that makes the leptonic top have a certain value of mass. Returns TRUE if there were real solutions. Returns FALSE if there were only complex solutions. In case of complex solutions, the real components of the solutions are given.

Parameters:
evInput:The event to solve.
tmassInput: The desired value of top quark mass in GeV.
nuz1Output: The first solution (smaller absolute value).
nuz2Output: The second solution.

Definition at line 82 of file Top_Decaykin.cc.

References a, abs, alpha, b, trackerHits::c, createTree::dd, hitfit::lepb_label, EgammaValidation_Wenu_cff::leptons, hitfit::Lepjets_Event::met(), mathSSE::sqrt(), hitfit::Lepjets_Event::sum(), and swap().

Referenced by hitfit::Top_Fit::fit_one_perm().

{
  bool discrim_flag = true;

  const Fourvec& vnu = ev.met ();
  Fourvec cprime = leptons (ev) + ev.sum (lepb_label);
  double alpha1 = tmass*tmass - cprime.m2();
  double a = 2 * 4 * (cprime.z()*cprime.z() - cprime.e()*cprime.e());
  double alpha = alpha1 + 2*(cprime.x()*vnu.x() + cprime.y()*vnu.y());
  double b = 4 * alpha * cprime.z();
  double c = alpha*alpha - 4 * cprime.e()*cprime.e() * vnu.vect().perp2();
  double d = b*b - 2*a*c;
  if (d < 0) {
    discrim_flag = false;
    d = 0;
  }

  double dd = sqrt (d);
  nuz1 = (-b + dd)/a;
  nuz2 = (-b - dd)/a;
  if (abs (nuz1) > abs (nuz2))
    swap (nuz1, nuz2);

  return discrim_flag;
}