CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullLepKinSolver.h
Go to the documentation of this file.
1 #ifndef TtFullLepKinSolver_h
2 #define TtFullLepKinSolver_h
3 
6 
7 #include "TLorentzVector.h"
8 #include "TMath.h"
9 
10 class TF2;
11 
12 /*
13  \class TtFullLepKinSolver TtFullLepKinSolver.h "TopQuarkAnalysis/TopKinFitter/interface/TtFullLepKinSolver.h"
14 
15  \brief Class to calculate solutions for neutrino momenta in dileptonic ttbar events and related probability weights
16 
17  Class to calculate solutions for neutrino momenta in dileptonic ttbar events and related probability weights.
18  A fourth-order polynomial in p_x(nu) is used with coefficients that are functions of the top-quark mass.
19  If physical (non-imaginary) solutions are found, the neutrino momenta are compared to the expected neutrino spectrum
20  (from simulation) to obtain a probability weight for each solution.
21  This class is based on a code by Jan Valenta.
22 
23 **/
24 
26 
27  public:
28 
31  double weight;
34  };
35 
39  TtFullLepKinSolver(const double, const double, const double, const std::vector<double>&, const double=80.4, const double=4.8);
42 
44  inline void useWeightFromMC(bool useMC) { useMCforBest_ = useMC; }
48  void SetConstraints(const double xx=0, const double yy=0);
50  NeutrinoSolution getNuSolution(const TLorentzVector& LV_l,
51  const TLorentzVector& LV_l_,
52  const TLorentzVector& LV_b,
53  const TLorentzVector& LV_b_);
54 
55  private:
56 
58  void FindCoeff(const TLorentzVector& al,
59  const TLorentzVector& l,
60  const TLorentzVector& b_al,
61  const TLorentzVector& b_l,
62  const double mt, const double mat, const double pxboost, const double pyboost,
63  double* q_coeff);
65  void TopRec(const TLorentzVector& al,
66  const TLorentzVector& l,
67  const TLorentzVector& b_al,
68  const TLorentzVector& b_l, const double sol);
70  double WeightSolfromMC() const;
72  double WeightSolfromShape() const;
74  int quartic(double* q_coeff, double* q_sol) const;
76  int cubic(const double* c_coeff, double* c_sol) const;
78  double sqr(const double x) const {return (x*x);}
80  void SWAP(double& realone, double& realtwo) const;
81 
82  private:
84  const double topmass_begin;
86  const double topmass_end;
88  const double topmass_step;
90  const double mw;
92  const double mb;
94  double pxmiss_, pymiss_;
95 
96  double C;
97  double D;
98  double F;
99  double pom;
100  double k16;
101  double k26;
102  double k36;
103  double k46;
104  double k56;
105  double k51;
106  double k61;
107  double m1;
108  double m2;
109  double m3;
110  double n1;
111  double n2;
112  double n3;
113 
115  TLorentzVector LV_n, LV_n_, LV_t, LV_t_, LV_tt_t, LV_tt_t_;
117  TLorentzVector genLV_n, genLV_n_;
118 
122  TF2* EventShape_;
123 };
124 
125 
126 #endif
int quartic(double *q_coeff, double *q_sol) const
int cubic(const double *c_coeff, double *c_sol) const
NeutrinoSolution getNuSolution(const TLorentzVector &LV_l, const TLorentzVector &LV_l_, const TLorentzVector &LV_b, const TLorentzVector &LV_b_)
double WeightSolfromMC() const
void SetConstraints(const double xx=0, const double yy=0)
const double topmass_end
void FindCoeff(const TLorentzVector &al, const TLorentzVector &l, const TLorentzVector &b_al, const TLorentzVector &b_l, const double mt, const double mat, const double pxboost, const double pyboost, double *q_coeff)
TtFullLepKinSolver()
default constructor
bool useMCforBest_
flag to swith from WeightSolfromMC() to WeightSolfromShape()
TLorentzVector genLV_n
provisional
const double topmass_step
void SWAP(double &realone, double &realtwo) const
TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution *asol)
TF2 * EventShape_
Event shape.
double sqr(const double x) const
TLorentzVector genLV_n_
void TopRec(const TLorentzVector &al, const TLorentzVector &l, const TLorentzVector &b_al, const TLorentzVector &b_l, const double sol)
void useWeightFromMC(bool useMC)
~TtFullLepKinSolver()
destructor
double WeightSolfromShape() const
use the parametrized event shape to obtain the solution weight.
TLorentzVector LV_tt_t
Definition: DDAxes.h:10
const double topmass_begin
TLorentzVector LV_tt_t_