CMS 3D CMS Logo

TtFullLepKinSolver Class Reference

#include <TopQuarkAnalysis/TopKinFitter/interface/TtFullLepKinSolver.h>

List of all members.

Public Member Functions

TtDilepEvtSolution addKinSolInfo (TtDilepEvtSolution *asol)
 TtFullLepKinSolver (double, double, double, double xx=0, double yy=0)
 TtFullLepKinSolver ()
void useWeightFromMC (bool useMC)
 ~TtFullLepKinSolver ()

Private Member Functions

int cubic (double *c_coeff, double *c_sol)
void FindCoeff (const TLorentzVector al, const TLorentzVector l, const TLorentzVector b_al, const TLorentzVector b_l, double mt, double mat, double pxboost, double pyboost, double *q_coeff)
int quartic (double *q_coeff, double *q_sol)
double sqr (double x)
void SWAP (double &realone, double &realtwo)
void TopRec (const TLorentzVector al, const TLorentzVector l, const TLorentzVector b_al, const TLorentzVector b_l, double sol)
double WeightSolfromMC ()
double WeightSolfromShape ()

Private Attributes

double C
double D
TF2 * EventShape_
double F
TLorentzVector genLV_n
TLorentzVector genLV_n_
double k16
double k26
double k36
double k46
double k51
double k56
double k61
TLorentzVector LV_n
TLorentzVector LV_n_
TLorentzVector LV_t
TLorentzVector LV_t_
TLorentzVector LV_tt_t
TLorentzVector LV_tt_t_
double m1
double m2
double m3
double mab
double maw
double mb
double mw
double n1
double n2
double n3
double pom
double pxmiss_
double pymiss_
double topmass_begin
double topmass_end
double topmass_step
bool useMCforBest_


Detailed Description

Definition at line 13 of file TtFullLepKinSolver.h.


Constructor & Destructor Documentation

TtFullLepKinSolver::TtFullLepKinSolver (  ) 

Definition at line 5 of file TtFullLepKinSolver.cc.

References EventShape_.

00006 {
00007   // That crude parametrisation has been obtained from a fit of O(1000) pythia events.
00008   // It is normalized to 1.
00009   EventShape_ = new TF2("landau2D","[0]*TMath::Landau(x,[1],[2],0)*TMath::Landau(y,[3],[4],0)",0,500,0,500);
00010   EventShape_->SetParameters(30.7137,56.2880,23.0744,59.1015,24.9145);
00011 }

TtFullLepKinSolver::TtFullLepKinSolver ( double  b,
double  e,
double  s,
double  xx = 0,
double  yy = 0 
)

Definition at line 13 of file TtFullLepKinSolver.cc.

References EventShape_, mab, maw, mb, mw, pxmiss_, pymiss_, topmass_begin, topmass_end, and topmass_step.

00014 {
00015   topmass_begin = b;
00016   topmass_end = e;
00017   topmass_step = s;
00018   pxmiss_ = xx;
00019   pymiss_ = yy;
00020   mw = 80.22;
00021   maw = 80.22;
00022   mab = 4.800;
00023   mb = 4.800;
00024   // That crude parametrisation has been obtained from a fit of O(1000) pythia events.
00025   // It is normalized to 1.
00026   EventShape_ = new TF2("landau2D","[0]*TMath::Landau(x,[1],[2],0)*TMath::Landau(y,[3],[4],0)",0,500,0,500);
00027   EventShape_->SetParameters(30.7137,56.2880,23.0744,59.1015,24.9145);
00028 }

TtFullLepKinSolver::~TtFullLepKinSolver (  ) 

Definition at line 33 of file TtFullLepKinSolver.cc.

References EventShape_.

00034 {
00035   delete EventShape_;
00036 }


Member Function Documentation

TtDilepEvtSolution TtFullLepKinSolver::addKinSolInfo ( TtDilepEvtSolution asol  ) 

Definition at line 38 of file TtFullLepKinSolver.cc.

References reco::Particle::energy(), FindCoeff(), genLV_n, genLV_n_, TtDilepEvtSolution::getCalJetB(), TtDilepEvtSolution::getCalJetBbar(), TtDilepEvtSolution::getElectronm(), TtDilepEvtSolution::getElectronp(), TtDilepEvtSolution::getGenN(), TtDilepEvtSolution::getGenNbar(), TtDilepEvtSolution::getMuonm(), TtDilepEvtSolution::getMuonp(), TtDilepEvtSolution::getTaum(), TtDilepEvtSolution::getTaup(), TtDilepEvtSolution::getWmDecay(), TtDilepEvtSolution::getWpDecay(), reco::Particle::px(), pxmiss_, reco::Particle::py(), pymiss_, reco::Particle::pz(), quartic(), TtDilepEvtSolution::setRecTopMass(), TtDilepEvtSolution::setRecWeightMax(), topmass_begin, topmass_end, topmass_step, TopRec(), useMCforBest_, weight, WeightSolfromMC(), and WeightSolfromShape().

Referenced by TtDilepEvtSolutionMaker::produce().

00039 {
00040   TtDilepEvtSolution fitsol(*asol);
00041   
00042   //antilepton and lepton
00043   TLorentzVector LV_e, LV_e_;
00044   //b and bbar quark
00045   TLorentzVector LV_b, LV_b_;
00046   
00047   bool hasMCinfo = true;
00048   if(fitsol.getGenN()) { // protect against non-dilept genevents
00049     genLV_n = TLorentzVector(fitsol.getGenN()->px(), fitsol.getGenN()->py(),
00050                              fitsol.getGenN()->pz(), fitsol.getGenN()->energy());
00051   } else hasMCinfo = false;
00052 
00053   if(fitsol.getGenNbar()) { // protect against non-dilept genevents
00054     genLV_n_ = TLorentzVector(fitsol.getGenNbar()->px(), fitsol.getGenNbar()->py(),
00055                               fitsol.getGenNbar()->pz(), fitsol.getGenNbar()->energy());
00056   } else hasMCinfo = false;
00057   // if MC is to be used to select the best top mass and is not available,
00058   // then nothing can be done. Stop here.
00059   if(useMCforBest_&&!hasMCinfo) return fitsol;
00060 
00061   // first lepton
00062   if (fitsol.getWpDecay() == "muon") {
00063     LV_e = TLorentzVector(fitsol.getMuonp().px(), fitsol.getMuonp().py(),
00064                           fitsol.getMuonp().pz(), fitsol.getMuonp().energy());
00065   } else if (fitsol.getWpDecay() == "electron") {
00066     LV_e = TLorentzVector(fitsol.getElectronp().px(), fitsol.getElectronp().py(),
00067                           fitsol.getElectronp().pz(), fitsol.getElectronp().energy());
00068   } else if (fitsol.getWpDecay() == "tau") {
00069     LV_e = TLorentzVector(fitsol.getTaup().px(), fitsol.getTaup().py(),
00070                           fitsol.getTaup().pz(), fitsol.getTaup().energy());
00071   }
00072     
00073   // second lepton
00074   if (fitsol.getWmDecay() == "muon") {
00075     LV_e_ = TLorentzVector(fitsol.getMuonm().px(), fitsol.getMuonm().py(),
00076                           fitsol.getMuonm().pz(), fitsol.getMuonm().energy());
00077   } else if (fitsol.getWmDecay() == "electron") {
00078     LV_e_ = TLorentzVector(fitsol.getElectronm().px(), fitsol.getElectronm().py(),
00079                            fitsol.getElectronm().pz(), fitsol.getElectronm().energy());
00080   } else if (fitsol.getWmDecay() == "tau") {
00081     LV_e_ = TLorentzVector(fitsol.getTaum().px(), fitsol.getTaum().py(),
00082                            fitsol.getTaum().pz(), fitsol.getTaum().energy());
00083   }
00084 
00085   // first jet
00086   LV_b = TLorentzVector(fitsol.getCalJetB().px(), fitsol.getCalJetB().py(),
00087                         fitsol.getCalJetB().pz(), fitsol.getCalJetB().energy());
00088 
00089   // second jet
00090   LV_b_ = TLorentzVector(fitsol.getCalJetBbar().px(), fitsol.getCalJetBbar().py(),
00091                          fitsol.getCalJetBbar().pz(), fitsol.getCalJetBbar().energy());
00092   
00093   //loop on top mass parameter
00094   double weightmax = -1e30;
00095   double mtmax = 0;
00096   for (double mt = topmass_begin; 
00097               mt < topmass_end + 0.5*topmass_step; 
00098               mt += topmass_step) {
00099     //cout << "mt = " << mt << endl;
00100     double q_coeff[5], q_sol[4];
00101     FindCoeff(LV_e, LV_e_, LV_b, LV_b_, mt, mt, pxmiss_, pymiss_, q_coeff);
00102     int NSol = quartic(q_coeff, q_sol);
00103     
00104     //loop on all solutions
00105     for (int isol = 0; isol < NSol; isol++) {
00106       TopRec(LV_e, LV_e_, LV_b, LV_b_, q_sol[isol]);
00107       double weight = useMCforBest_ ? WeightSolfromMC() : WeightSolfromShape();
00108       if (weight > weightmax) {
00109         weightmax =weight;
00110         mtmax = mt;
00111       }
00112     }
00113     
00114     //for (int i=0;i<5;i++) cout << " q_coeff["<<i<< "]= " << q_coeff[i];
00115     //cout << endl;
00116     
00117     //for (int i=0;i<4;i++) cout << " q_sol["<<i<< "]= " << q_sol[i];
00118     //cout << endl;
00119     //cout << "NSol_" << NSol << endl;
00120   }
00121   
00122   fitsol.setRecTopMass(mtmax);
00123   fitsol.setRecWeightMax(weightmax);
00124   
00125   return fitsol;
00126 }

int TtFullLepKinSolver::cubic ( double *  c_coeff,
double *  c_sol 
) [private]

Definition at line 331 of file TtFullLepKinSolver.cc.

References funct::cos(), h, i, p, phi, Pi, sqr(), funct::sqrt(), SWAP(), and w.

Referenced by quartic().

00331                                                             {
00332   int i, nreal;
00333   double w, p, q, dis, h, phi;
00334   
00335   if (coeffs[3]!=0.0) {
00336     /* cubic problem? */
00337     w = coeffs[2]/(3*coeffs[3]);
00338     p = sqr(coeffs[1]/(3*coeffs[3])-sqr(w))*(coeffs[1]/(3*coeffs[3])-sqr(w));
00339     q = -0.5*(2*sqr(w)*w-(coeffs[1]*w-coeffs[0])/coeffs[3]);
00340     dis = sqr(q)+p;
00341     /* discriminant */
00342     if (dis<0.0) {
00343       /* 3 real solutions */
00344       h = q/sqrt(-p);
00345       if (h>1.0) h = 1.0;
00346       /* confine the argument of */
00347       if (h<-1.0) h = -1.0;
00348       /* acos to [-1;+1] */
00349       phi = acos(h);
00350       p = 2*TMath::Power(-p, 1.0/6.0);
00351       for (i=0; i<3; i++) 
00352         koreny[i] = p*cos((phi+2*i*TMath::Pi())/3.0) - w;
00353       if (koreny[1]<koreny[0]) SWAP(koreny[0], koreny[1]);
00354       /* sort results */
00355       if (koreny[2]<koreny[1]) SWAP(koreny[1], koreny[2]);
00356       if (koreny[1]<koreny[0]) SWAP(koreny[0], koreny[1]);
00357       nreal = 3;
00358     }
00359     else {
00360       /* only one real solution */
00361       dis = sqrt(dis);
00362       h = TMath::Power(fabs(q+dis), 1.0/3.0);
00363       p = TMath::Power(fabs(q-dis), 1.0/3.0);
00364       koreny[0] = ((q+dis>0.0)? h : -h) + ((q-dis>0.0)? p : -p) -  w;
00365       nreal = 1;
00366     }
00367 
00368     /* Perform one step of a Newton iteration in order to minimize
00369        round-off errors */
00370     for (i=0; i<nreal; i++) {
00371       h = coeffs[1] + koreny[i] * (2 * coeffs[2] + 3 * koreny[i] * coeffs[3]);
00372       if (h != 0.0)
00373         koreny[i] -= (coeffs[0] + koreny[i] * (coeffs[1] + koreny[i] * (coeffs[2] + koreny[i] * coeffs[3])))/h;
00374     }
00375   }
00376 
00377   else if (coeffs[2]!=0.0) {
00378     /* quadratic problem? */
00379     p = 0.5*coeffs[1]/coeffs[2];
00380     dis = sqr(p) - coeffs[0]/coeffs[2];
00381     if (dis>=0.0) {
00382       /* two real solutions */
00383       dis = sqrt(dis);
00384       koreny[0] = -p - dis;
00385       koreny[1] = -p + dis;
00386       nreal = 2;
00387     }
00388     else
00389       /* no real solution */
00390       nreal = 0;
00391   }
00392 
00393   else if (coeffs[1]!=0.0) {
00394     /* linear problem? */
00395     koreny[0] = -coeffs[0]/coeffs[1];
00396     nreal = 1;
00397   }
00398 
00399   else
00400     /* no equation */
00401     nreal = 0;
00402   
00403   return nreal;
00404 }

void TtFullLepKinSolver::FindCoeff ( const TLorentzVector  al,
const TLorentzVector  l,
const TLorentzVector  b_al,
const TLorentzVector  b_l,
double  mt,
double  mat,
double  pxboost,
double  pyboost,
double *  q_coeff 
) [private]

Definition at line 128 of file TtFullLepKinSolver.cc.

References C, D, F, int, k16, k26, k36, k46, k51, k56, k61, l1, l2, l3, m1, m2, m3, mab, maw, mb, mw, n1, n2, n3, pom, funct::pow(), and sqr().

Referenced by addKinSolInfo().

00133                                                     {
00134 
00135 double E, apom1, apom2, apom3;
00136 double k11, k21, k31, k41,  cpom1, cpom2, cpom3, l11, l21, l31, l41, l51, l61, k1, k2, k3, k4, k5,k6;
00137 double l1, l2, l3, l4, l5, l6, k15, k25, k35, k45;
00138 
00139   C = -al.Px()-b_al.Px()-l.Px()-b_l.Px() + px_miss;
00140   D = -al.Py()-b_al.Py()-l.Py()-b_l.Py() + py_miss;
00141 
00142   // right side offirst two linear equations - missing pT
00143 
00144   
00145   E = (sqr(mt)-sqr(mw)-sqr(mb))/(2*b_al.E())-sqr(mw)/(2*al.E())-al.E()+al.Px()*b_al.Px()/b_al.E()+al.Py()*b_al.Py()/b_al.E()+al.Pz()*b_al.Pz()/b_al.E();
00146   F = (sqr(mat)-sqr(maw)-sqr(mab))/(2*b_l.E())-sqr(maw)/(2*l.E())-l.E()+l.Px()*b_l.Px()/b_l.E()+l.Py()*b_l.Py()/b_l.E()+l.Pz()*b_l.Pz()/b_l.E();
00147   
00148   m1 = al.Px()/al.E()-b_al.Px()/b_al.E();
00149   m2 = al.Py()/al.E()-b_al.Py()/b_al.E();
00150   m3 = al.Pz()/al.E()-b_al.Pz()/b_al.E();
00151   
00152   n1 = l.Px()/l.E()-b_l.Px()/b_l.E();
00153   n2 = l.Py()/l.E()-b_l.Py()/b_l.E();
00154   n3 = l.Pz()/l.E()-b_l.Pz()/b_l.E();
00155   
00156   pom = E-m1*C-m2*D;
00157   apom1 = sqr(al.Px())-sqr(al.E());
00158   apom2 = sqr(al.Py())-sqr(al.E());
00159   apom3 = sqr(al.Pz())-sqr(al.E());
00160   
00161   k11 = 1/sqr(al.E())*(pow(mw,4)/4+sqr(C)*apom1+sqr(D)*apom2+apom3*sqr(pom)/sqr(m3)+sqr(mw)*(al.Px()*C+al.Py()*D+al.Pz()*pom/m3)+2*al.Px()*al.Py()*C*D+2*al.Px()*al.Pz()*C*pom/m3+2*al.Py()*al.Pz()*D*pom/m3);
00162   k21 = 1/sqr(al.E())*(-2*C*m3*n3*apom1+2*apom3*n3*m1*pom/m3-sqr(mw)*m3*n3*al.Px()+sqr(mw)*m1*n3*al.Pz()-2*al.Px()*al.Py()*D*m3*n3+2*al.Px()*al.Pz()*C*m1*n3-2*al.Px()*al.Pz()*n3*pom+2*al.Py()*al.Pz()*D*m1*n3);
00163   k31 = 1/sqr(al.E())*(-2*D*m3*n3*apom2+2*apom3*n3*m2*pom/m3-sqr(mw)*m3*n3*al.Py()+sqr(mw)*m2*n3*al.Pz()-2*al.Px()*al.Py()*C*m3*n3+2*al.Px()*al.Pz()*C*m2*n3-2*al.Py()*al.Pz()*n3*pom+2*al.Py()*al.Pz()*D*m2*n3);
00164   k41 = 1/sqr(al.E())*(2*apom3*m1*m2*sqr(n3)+2*al.Px()*al.Py()*sqr(m3)*sqr(n3)-2*al.Px()*al.Pz()*m2*m3*sqr(n3)-2*al.Py()*al.Pz()*m1*m3*sqr(n3));
00165   k51 = 1/sqr(al.E())*(apom1*sqr(m3)*sqr(n3)+apom3*sqr(m1)*sqr(n3)-2*al.Px()*al.Pz()*m1*m3*sqr(n3));
00166   k61 = 1/sqr(al.E())*(apom2*sqr(m3)*sqr(n3)+apom3*sqr(m2)*sqr(n3)-2*al.Py()*al.Pz()*m2*m3*sqr(n3));
00167   
00168   cpom1 = sqr(l.Px())-sqr(l.E());
00169   cpom2 = sqr(l.Py())-sqr(l.E());
00170   cpom3 = sqr(l.Pz())-sqr(l.E());
00171   
00172   l11 = 1/sqr(l.E())*(pow(maw,4)/4+cpom3*sqr(F)/sqr(n3)+sqr(maw)*l.Pz()*F/n3);
00173   l21 = 1/sqr(l.E())*(-2*cpom3*F*m3*n1/n3+sqr(maw)*(l.Px()*m3*n3-l.Pz()*n1*m3)+2*l.Px()*l.Pz()*F*m3);
00174   l31 = 1/sqr(l.E())*(-2*cpom3*F*m3*n2/n3+sqr(maw)*(l.Py()*m3*n3-l.Pz()*n2*m3)+2*l.Py()*l.Pz()*F*m3);
00175   l41 = 1/sqr(l.E())*(2*cpom3*n1*n2*sqr(m3)+2*l.Px()*l.Py()*sqr(m3)*sqr(n3)-2*l.Px()*l.Pz()*n2*n3*sqr(m3)-2*l.Py()*l.Pz()*n1*n3*sqr(m3));
00176   l51 = 1/sqr(l.E())*(cpom1*sqr(m3)*sqr(n3)+cpom3*sqr(n1)*sqr(m3)-2*l.Px()*l.Pz()*n1*n3*sqr(m3));
00177   l61 = 1/sqr(l.E())*(cpom2*sqr(m3)*sqr(n3)+cpom3*sqr(n2)*sqr(m3)-2*l.Py()*l.Pz()*n2*n3*sqr(m3));
00178   
00179   k1 = k11*k61;
00180   k2 = k61*k21/k51;
00181   k3 = k31;
00182   k4 = k41/k51;
00183   k5 = k61/k51;
00184   k6 = 1;
00185   
00186   l1 = l11*k61;
00187   l2 = l21*k61/k51;
00188   l3 = l31;
00189   l4 = l41/k51;
00190   l5 = l51*k61/(sqr(k51));
00191   l6 = l61/k61;
00192   
00193   k15 = k1*l5-l1*k5;
00194   k25 = k2*l5-l2*k5;
00195   k35 = k3*l5-l3*k5;
00196   k45 = k4*l5-l4*k5;
00197   
00198   k16 = k1*l6-l1*k6;
00199   k26 = k2*l6-l2*k6;
00200   k36 = k3*l6-l3*k6;
00201   k46 = k4*l6-l4*k6;
00202   k56 = k5*l6-l5*k6;
00203   
00204 
00205   koeficienty[0] = k15*sqr(k36)-k35*k36*k16-k56*sqr(k16);
00206   koeficienty[1] = 2*k15*k36*k46+k25*sqr(k36)+k35*(-k46*k16-k36*k26)-k45*k36*k16-2*k56*k26*k16;
00207   koeficienty[2] = k15*sqr(k46)+2*k25*k36*k46+k35*(-k46*k26-k36*k56)-k56*(sqr(k26)+2*k56*k16)-k45*(k46*k16+k36*k26);
00208   koeficienty[3] = k25*sqr(k46)-k35*k46*k56-k45*(k46*k26+k36*k56)-2*sqr(k56)*k26;
00209   koeficienty[4] = -k45*k46*k56-pow(k56,3);
00210   
00211   // normalization of coefficients
00212   int moc=(int(log10(fabs(koeficienty[0])))+int(log10(fabs(koeficienty[4]))))/2;
00213   
00214   koeficienty[0]=koeficienty[0]/TMath::Power(10,moc);
00215   koeficienty[1]=koeficienty[1]/TMath::Power(10,moc);
00216   koeficienty[2]=koeficienty[2]/TMath::Power(10,moc);
00217   koeficienty[3]=koeficienty[3]/TMath::Power(10,moc);
00218   koeficienty[4]=koeficienty[4]/TMath::Power(10,moc);
00219 }

int TtFullLepKinSolver::quartic ( double *  q_coeff,
double *  q_sol 
) [private]

Definition at line 270 of file TtFullLepKinSolver.cc.

References b1, b2, c, cubic(), debug_cff::d0, d1, h, i, sqr(), funct::sqrt(), t, w, and z.

Referenced by addKinSolInfo().

00270                                                                    {
00271   int i, nreal;
00272   double w, b0, b1, b2;
00273   double c[4];
00274   double d0, d1, h, t, z;
00275   double *px;
00276 
00277   
00278   if (koeficienty[4]==0.0) 
00279     return cubic(koeficienty, koreny);
00280   /* quartic problem? */
00281   w = koeficienty[3]/(4*koeficienty[4]);
00282   /* offset */
00283   b2 = -6*sqr(w) + koeficienty[2]/koeficienty[4];
00284   /* koeficienty. of shifted polynomial */
00285   b1 = (8*sqr(w) - 2*koeficienty[2]/koeficienty[4])*w + koeficienty[1]/koeficienty[4];
00286   b0 = ((-3*sqr(w) + koeficienty[2]/koeficienty[4])*w - koeficienty[1]/koeficienty[4])*w + koeficienty[0]/koeficienty[4];
00287 
00288   c[3] = 1.0;
00289   /* cubic resolvent */
00290   c[2] = b2;
00291   c[1] = -4*b0;
00292   c[0] = sqr(b1) - 4*b0*b2;
00293   
00294   
00295   i = cubic(c, koreny);
00296   z = koreny[0];
00297   //double z1=1.0,z2=2.0,z3=3.0;
00298   //TMath::RootsCubic(c,z1,z2,z3);
00299   //if (z2 !=0) z = z2;
00300   //if (z1 !=0) z = z1;
00301   /* only lowermost root needed */
00302 
00303   nreal = 0;
00304   px = koreny;
00305   t = sqrt(0.25*sqr(z) - b0);
00306   for (i=-1; i<=1; i+=2) {
00307     d0 = -0.5*z + i*t;
00308     /* coeffs. of quadratic factor */
00309     d1 = (t!=0.0)? -i*0.5*b1/t : i*sqrt(-z - b2);
00310     h = 0.25*sqr(d1) - d0;
00311     if (h>=0.0) {
00312       h = sqrt(h);
00313       nreal += 2;
00314       *px++ = -0.5*d1 - h - w;
00315       *px++ = -0.5*d1 + h - w;
00316     }
00317   }
00318 
00319   if (nreal==4) {
00320     /* sort results */
00321 //    if (koreny[2]<koreny[0]) SWAP(koreny[0], koreny[2]);
00322 //    if (koreny[3]<koreny[1]) SWAP(koreny[1], koreny[3]);
00323 //    if (koreny[1]<koreny[0]) SWAP(koreny[0], koreny[1]);
00324 //    if (koreny[3]<koreny[2]) SWAP(koreny[2], koreny[3]);
00325 //    if (koreny[2]<koreny[1]) SWAP(koreny[1], koreny[2]);
00326   }
00327   return nreal;
00328 
00329 }

double TtFullLepKinSolver::sqr ( double  x  )  [inline, private]

Definition at line 41 of file TtFullLepKinSolver.h.

Referenced by cubic(), FindCoeff(), and quartic().

00041 {return (x*x);}

void TtFullLepKinSolver::SWAP ( double &  realone,
double &  realtwo 
) [private]

Definition at line 407 of file TtFullLepKinSolver.cc.

References edm::aux.

Referenced by cubic().

00407                                                                {
00408   if (realtwo < realone) {
00409     double aux = realtwo;
00410     realtwo = realone;
00411     realone = aux;
00412   }
00413 }

void TtFullLepKinSolver::TopRec ( const TLorentzVector  al,
const TLorentzVector  l,
const TLorentzVector  b_al,
const TLorentzVector  b_l,
double  sol 
) [private]

Definition at line 221 of file TtFullLepKinSolver.cc.

References edm::aux, C, D, F, k16, k26, k36, k46, k51, k56, k61, LV_n, LV_n_, LV_t, LV_t_, LV_tt_t, LV_tt_t_, m1, m2, m3, n1, n2, n3, pom, funct::pow(), and muonGeometry::TVector3().

Referenced by addKinSolInfo().

00224                                                                   {
00225 
00226   TVector3 t_ttboost;
00227   TLorentzVector aux;
00228   double pxp, pyp, pzp, pup, pvp, pwp;
00229   
00230     
00231   pxp = sol*(m3*n3/k51);   
00232   pyp = -(m3*n3/k61)*(k56*pow(sol,2) + k26*sol + k16)/(k36 + k46*sol);
00233   pzp = -1/n3*(n1*pxp + n2*pyp - F);
00234   pwp = 1/m3*(m1*pxp + m2*pyp + pom);
00235   pup = C - pxp;
00236   pvp = D - pyp;
00237      
00238   LV_n_.SetXYZM(pxp, pyp, pzp, 0.0);
00239   LV_n.SetXYZM(pup, pvp, pwp, 0.0);
00240       
00241   
00242   LV_t_ = b_l + l + LV_n_;
00243   LV_t = b_al + al + LV_n;  
00244   
00245   
00246   aux = (LV_t_ + LV_t);
00247   t_ttboost = -aux.BoostVector();
00248   LV_tt_t_ = LV_t_;
00249   LV_tt_t = LV_t;
00250   LV_tt_t_.Boost(t_ttboost);
00251   LV_tt_t.Boost(t_ttboost); 
00252 }

void TtFullLepKinSolver::useWeightFromMC ( bool  useMC  )  [inline]

Definition at line 21 of file TtFullLepKinSolver.h.

References useMCforBest_.

Referenced by TtDilepEvtSolutionMaker::produce().

00021 { useMCforBest_ = useMC; }

double TtFullLepKinSolver::WeightSolfromMC (  )  [private]

Definition at line 254 of file TtFullLepKinSolver.cc.

References genLV_n, genLV_n_, LV_n, LV_n_, and weight.

Referenced by addKinSolInfo().

00254                                            {
00255 
00256   double weight = 1;
00257   weight = ((LV_n.E() > genLV_n.E())? genLV_n.E()/LV_n.E(): LV_n.E()/genLV_n.E())
00258            *((LV_n_.E() > genLV_n_.E())? genLV_n_.E()/LV_n_.E(): LV_n_.E()/genLV_n_.E());
00259   return weight;
00260 
00261 }

double TtFullLepKinSolver::WeightSolfromShape (  )  [private]

Definition at line 263 of file TtFullLepKinSolver.cc.

References EventShape_, LV_n, and LV_n_.

Referenced by addKinSolInfo().

00263                                               {
00264   
00265   // Use the parametrized event shape to obtain the solution weight.
00266   return EventShape_->Eval(LV_n.E(),LV_n_.E());
00267 
00268 }


Member Data Documentation

double TtFullLepKinSolver::C [private]

Definition at line 52 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::D [private]

Definition at line 53 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

TF2* TtFullLepKinSolver::EventShape_ [private]

Definition at line 77 of file TtFullLepKinSolver.h.

Referenced by TtFullLepKinSolver(), WeightSolfromShape(), and ~TtFullLepKinSolver().

double TtFullLepKinSolver::F [private]

Definition at line 54 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

TLorentzVector TtFullLepKinSolver::genLV_n [private]

Definition at line 73 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and WeightSolfromMC().

TLorentzVector TtFullLepKinSolver::genLV_n_ [private]

Definition at line 73 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and WeightSolfromMC().

double TtFullLepKinSolver::k16 [private]

Definition at line 56 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::k26 [private]

Definition at line 57 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::k36 [private]

Definition at line 58 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::k46 [private]

Definition at line 59 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::k51 [private]

Definition at line 61 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::k56 [private]

Definition at line 60 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::k61 [private]

Definition at line 62 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

TLorentzVector TtFullLepKinSolver::LV_n [private]

Definition at line 70 of file TtFullLepKinSolver.h.

Referenced by TopRec(), WeightSolfromMC(), and WeightSolfromShape().

TLorentzVector TtFullLepKinSolver::LV_n_ [private]

Definition at line 70 of file TtFullLepKinSolver.h.

Referenced by TopRec(), WeightSolfromMC(), and WeightSolfromShape().

TLorentzVector TtFullLepKinSolver::LV_t [private]

Definition at line 70 of file TtFullLepKinSolver.h.

Referenced by TopRec().

TLorentzVector TtFullLepKinSolver::LV_t_ [private]

Definition at line 70 of file TtFullLepKinSolver.h.

Referenced by TopRec().

TLorentzVector TtFullLepKinSolver::LV_tt_t [private]

Definition at line 70 of file TtFullLepKinSolver.h.

Referenced by TopRec().

TLorentzVector TtFullLepKinSolver::LV_tt_t_ [private]

Definition at line 70 of file TtFullLepKinSolver.h.

Referenced by TopRec().

double TtFullLepKinSolver::m1 [private]

Definition at line 63 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::m2 [private]

Definition at line 64 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::m3 [private]

Definition at line 65 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::mab [private]

Definition at line 50 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TtFullLepKinSolver().

double TtFullLepKinSolver::maw [private]

Definition at line 50 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TtFullLepKinSolver().

double TtFullLepKinSolver::mb [private]

Definition at line 50 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TtFullLepKinSolver().

double TtFullLepKinSolver::mw [private]

Definition at line 50 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TtFullLepKinSolver().

double TtFullLepKinSolver::n1 [private]

Definition at line 66 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::n2 [private]

Definition at line 67 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::n3 [private]

Definition at line 68 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::pom [private]

Definition at line 55 of file TtFullLepKinSolver.h.

Referenced by FindCoeff(), and TopRec().

double TtFullLepKinSolver::pxmiss_ [private]

Definition at line 48 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and TtFullLepKinSolver().

double TtFullLepKinSolver::pymiss_ [private]

Definition at line 48 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and TtFullLepKinSolver().

double TtFullLepKinSolver::topmass_begin [private]

Definition at line 45 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and TtFullLepKinSolver().

double TtFullLepKinSolver::topmass_end [private]

Definition at line 46 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and TtFullLepKinSolver().

double TtFullLepKinSolver::topmass_step [private]

Definition at line 47 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and TtFullLepKinSolver().

bool TtFullLepKinSolver::useMCforBest_ [private]

Definition at line 75 of file TtFullLepKinSolver.h.

Referenced by addKinSolInfo(), and useWeightFromMC().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:34:43 2009 for CMSSW by  doxygen 1.5.4