CMS 3D CMS Logo

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 public:
29  double weight;
32  };
33 
38  const double, const double, const double, const std::vector<double>&, const double = 80.4, const double = 4.8);
41 
43  inline void useWeightFromMC(bool useMC) { useMCforBest_ = useMC; }
47  void SetConstraints(const double xx = 0, const double yy = 0);
49  NeutrinoSolution getNuSolution(const TLorentzVector& LV_l,
50  const TLorentzVector& LV_l_,
51  const TLorentzVector& LV_b,
52  const TLorentzVector& LV_b_);
53 
54 private:
56  void FindCoeff(const TLorentzVector& al,
57  const TLorentzVector& l,
58  const TLorentzVector& b_al,
59  const TLorentzVector& b_l,
60  const double mt,
61  const double mat,
62  const double pxboost,
63  const double pyboost,
64  double* q_coeff);
66  void TopRec(const TLorentzVector& al,
67  const TLorentzVector& l,
68  const TLorentzVector& b_al,
69  const TLorentzVector& b_l,
70  const double sol);
72  double WeightSolfromMC() const;
74  double WeightSolfromShape() const;
76  int quartic(double* q_coeff, double* q_sol) const;
78  int cubic(const double* c_coeff, double* c_sol) const;
80  double sqr(const double x) const { return (x * x); }
82  void SWAP(double& realone, double& realtwo) const;
83 
84 private:
86  const double topmass_begin;
88  const double topmass_end;
90  const double topmass_step;
92  const double mw;
94  const double mb;
96  double pxmiss_, pymiss_;
97 
98  double C;
99  double D;
100  double F;
101  double pom;
102  double k16;
103  double k26;
104  double k36;
105  double k46;
106  double k56;
107  double k51;
108  double k61;
109  double m1;
110  double m2;
111  double m3;
112  double n1;
113  double n2;
114  double n3;
115 
117  TLorentzVector LV_n, LV_n_, LV_t, LV_t_, LV_tt_t, LV_tt_t_;
119  TLorentzVector genLV_n, genLV_n_;
120 
125 };
126 
127 #endif
double WeightSolfromMC() const
NeutrinoSolution getNuSolution(const TLorentzVector &LV_l, const TLorentzVector &LV_l_, const TLorentzVector &LV_b, const TLorentzVector &LV_b_)
int cubic(const double *c_coeff, double *c_sol) 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()
double WeightSolfromShape() const
use the parametrized event shape to obtain the solution weight.
TLorentzVector genLV_n
provisional
const double topmass_step
void SWAP(double &realone, double &realtwo) const
dictionary useMC
Definition: JetHT_cfg.py:45
TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution *asol)
int quartic(double *q_coeff, double *q_sol) const
TF2 * EventShape_
Event shape.
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
constexpr float sol
Definition: Config.h:13
double sqr(const double x) const
TLorentzVector LV_tt_t
const double topmass_begin
TLorentzVector LV_tt_t_