CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/TopQuarkAnalysis/TopKinFitter/interface/TtFullLepKinSolver.h

Go to the documentation of this file.
00001 #ifndef TtFullLepKinSolver_h
00002 #define TtFullLepKinSolver_h
00003 
00004 #include "AnalysisDataFormats/TopObjects/interface/TtDilepEvtSolution.h"
00005 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00006 
00007 #include "TLorentzVector.h"
00008 #include "TMath.h"
00009 
00010 class TF2;
00011 
00012 /*
00013   \class   TtFullLepKinSolver TtFullLepKinSolver.h "TopQuarkAnalysis/TopKinFitter/interface/TtFullLepKinSolver.h"
00014   
00015   \brief   Class to calculate solutions for neutrino momenta in dileptonic ttbar events and related probability weights
00016 
00017   Class to calculate solutions for neutrino momenta in dileptonic ttbar events and related probability weights.
00018   A fourth-order polynomial in p_x(nu) is used with coefficients that are functions of the top-quark mass.
00019   If physical (non-imaginary) solutions are found, the neutrino momenta are compared to the expected neutrino spectrum
00020   (from simulation) to obtain a probability weight for each solution.
00021   This class is based on a code by Jan Valenta.
00022   
00023 **/
00024 
00025 class TtFullLepKinSolver {
00026 
00027  public:
00028 
00030   struct NeutrinoSolution {
00031     double weight;
00032     reco::LeafCandidate neutrino;
00033     reco::LeafCandidate neutrinoBar; 
00034   };
00035 
00037   TtFullLepKinSolver();
00039   TtFullLepKinSolver(const double, const double, const double, const std::vector<double>, const double=80.4, const double=4.8);
00041   ~TtFullLepKinSolver();
00042 
00044   inline void useWeightFromMC(bool useMC) { useMCforBest_ = useMC; }
00046   TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution * asol); 
00048   void SetConstraints(const double xx=0, const double yy=0);
00050   NeutrinoSolution getNuSolution(TLorentzVector LV_l, 
00051                                  TLorentzVector LV_l_, 
00052                                  TLorentzVector LV_b, 
00053                                  TLorentzVector LV_b_);                           
00054                              
00055  private:
00056 
00058   void FindCoeff(const TLorentzVector al, 
00059                  const TLorentzVector l,
00060                  const TLorentzVector b_al,
00061                  const TLorentzVector b_l,
00062                  const double mt, const double mat, const double pxboost, const double pyboost,
00063                  double* q_coeff);
00065   void TopRec(const TLorentzVector al, 
00066               const TLorentzVector l,
00067               const TLorentzVector b_al,
00068               const TLorentzVector b_l, const double sol);
00070   double WeightSolfromMC() const;
00072   double WeightSolfromShape() const;
00074   int quartic(double* q_coeff, double* q_sol) const;
00076   int cubic(const double* c_coeff, double* c_sol) const;
00078   double sqr(const double x) const {return (x*x);}
00080   void SWAP(double& realone, double& realtwo) const;
00081     
00082  private:
00084   const double topmass_begin;
00086   const double topmass_end;
00088   const double topmass_step;
00090   const double mw;
00092   const double mb;
00094   double pxmiss_, pymiss_;
00095   
00096   double C;
00097   double D;
00098   double F;
00099   double pom;
00100   double k16;
00101   double k26;
00102   double k36;
00103   double k46;
00104   double k56;
00105   double k51;
00106   double k61;
00107   double m1;
00108   double m2;
00109   double m3;
00110   double n1;
00111   double n2;
00112   double n3;
00113   
00115   TLorentzVector LV_n, LV_n_, LV_t, LV_t_, LV_tt_t, LV_tt_t_;  
00117   TLorentzVector genLV_n, genLV_n_;  
00118     
00120   bool useMCforBest_;
00122   TF2* EventShape_;  
00123 };
00124 
00125 
00126 #endif