CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullLepKinSolver.cc
Go to the documentation of this file.
2 #include "TF2.h"
3 
5  topmass_begin(0),
6  topmass_end(0),
7  topmass_step(0),
8  mw(80.4),
9  mb(4.8),
10  pxmiss_(0),
11  pymiss_(0)
12 {
13  // That crude parametrisation has been obtained from a fit of O(1000) pythia events.
14  // It is normalized to 1.
15  EventShape_ = new TF2("landau2D","[0]*TMath::Landau(x,[1],[2],0)*TMath::Landau(y,[3],[4],0)",0,500,0,500);
16  EventShape_->SetParameters(30.7137,56.2880,23.0744,59.1015,24.9145);
17 }
18 
19 TtFullLepKinSolver::TtFullLepKinSolver(const double b, const double e, const double s, const std::vector<double>& nupars, const double mW, const double mB):
20  topmass_begin(b),
21  topmass_end(e),
22  topmass_step(s),
23  mw(mW),
24  mb(mB),
25  pxmiss_(0),
26  pymiss_(0)
27 {
28  EventShape_ = new TF2("landau2D","[0]*TMath::Landau(x,[1],[2],0)*TMath::Landau(y,[3],[4],0)",0,500,0,500);
29  EventShape_->SetParameters(nupars[0],nupars[1],nupars[2],nupars[3],nupars[4]);
30 }
31 
32 //
33 // destructor
34 //
36 {
37  delete EventShape_;
38 }
39 
41 {
42  TtDilepEvtSolution fitsol(*asol);
43 
44  //antilepton and lepton
45  TLorentzVector LV_e, LV_e_;
46  //b and bbar quark
47  TLorentzVector LV_b, LV_b_;
48 
49  bool hasMCinfo = true;
50  if(fitsol.getGenN()) { // protect against non-dilept genevents
51  genLV_n = TLorentzVector(fitsol.getGenN()->px(), fitsol.getGenN()->py(),
52  fitsol.getGenN()->pz(), fitsol.getGenN()->energy());
53  } else hasMCinfo = false;
54 
55  if(fitsol.getGenNbar()) { // protect against non-dilept genevents
56  genLV_n_ = TLorentzVector(fitsol.getGenNbar()->px(), fitsol.getGenNbar()->py(),
57  fitsol.getGenNbar()->pz(), fitsol.getGenNbar()->energy());
58  } else hasMCinfo = false;
59  // if MC is to be used to select the best top mass and is not available,
60  // then nothing can be done. Stop here.
61  if(useMCforBest_&&!hasMCinfo) return fitsol;
62 
63  // first lepton
64  if (fitsol.getWpDecay() == "muon") {
65  LV_e = TLorentzVector(fitsol.getMuonp().px(), fitsol.getMuonp().py(),
66  fitsol.getMuonp().pz(), fitsol.getMuonp().energy());
67  } else if (fitsol.getWpDecay() == "electron") {
68  LV_e = TLorentzVector(fitsol.getElectronp().px(), fitsol.getElectronp().py(),
69  fitsol.getElectronp().pz(), fitsol.getElectronp().energy());
70  } else if (fitsol.getWpDecay() == "tau") {
71  LV_e = TLorentzVector(fitsol.getTaup().px(), fitsol.getTaup().py(),
72  fitsol.getTaup().pz(), fitsol.getTaup().energy());
73  }
74 
75  // second lepton
76  if (fitsol.getWmDecay() == "muon") {
77  LV_e_ = TLorentzVector(fitsol.getMuonm().px(), fitsol.getMuonm().py(),
78  fitsol.getMuonm().pz(), fitsol.getMuonm().energy());
79  } else if (fitsol.getWmDecay() == "electron") {
80  LV_e_ = TLorentzVector(fitsol.getElectronm().px(), fitsol.getElectronm().py(),
81  fitsol.getElectronm().pz(), fitsol.getElectronm().energy());
82  } else if (fitsol.getWmDecay() == "tau") {
83  LV_e_ = TLorentzVector(fitsol.getTaum().px(), fitsol.getTaum().py(),
84  fitsol.getTaum().pz(), fitsol.getTaum().energy());
85  }
86 
87  // first jet
88  LV_b = TLorentzVector(fitsol.getCalJetB().px(), fitsol.getCalJetB().py(),
89  fitsol.getCalJetB().pz(), fitsol.getCalJetB().energy());
90 
91  // second jet
92  LV_b_ = TLorentzVector(fitsol.getCalJetBbar().px(), fitsol.getCalJetBbar().py(),
93  fitsol.getCalJetBbar().pz(), fitsol.getCalJetBbar().energy());
94 
95  //loop on top mass parameter
96  double weightmax = -1e30;
97  double mtmax = 0;
98  for (double mt = topmass_begin;
99  mt < topmass_end + 0.5*topmass_step;
100  mt += topmass_step) {
101  //cout << "mt = " << mt << endl;
102  double q_coeff[5], q_sol[4];
103  FindCoeff(LV_e, LV_e_, LV_b, LV_b_, mt, mt, pxmiss_, pymiss_, q_coeff);
104  int NSol = quartic(q_coeff, q_sol);
105 
106  //loop on all solutions
107  for (int isol = 0; isol < NSol; isol++) {
108  TopRec(LV_e, LV_e_, LV_b, LV_b_, q_sol[isol]);
110  if (weight > weightmax) {
111  weightmax =weight;
112  mtmax = mt;
113  }
114  }
115 
116  //for (int i=0;i<5;i++) cout << " q_coeff["<<i<< "]= " << q_coeff[i];
117  //cout << endl;
118 
119  //for (int i=0;i<4;i++) cout << " q_sol["<<i<< "]= " << q_sol[i];
120  //cout << endl;
121  //cout << "NSol_" << NSol << endl;
122  }
123 
124  fitsol.setRecTopMass(mtmax);
125  fitsol.setRecWeightMax(weightmax);
126 
127  return fitsol;
128 }
129 
130 void
131 TtFullLepKinSolver::SetConstraints(const double xx, const double yy)
132 {
133  pxmiss_ = xx;
134  pymiss_ = yy;
135 }
136 
138 TtFullLepKinSolver::getNuSolution(const TLorentzVector& LV_l,
139  const TLorentzVector& LV_l_,
140  const TLorentzVector& LV_b,
141  const TLorentzVector& LV_b_)
142 {
145 
146  //loop on top mass parameter
147  double weightmax = -1;
148  for(double mt = topmass_begin;
149  mt < topmass_end + 0.5*topmass_step;
150  mt += topmass_step) {
151  double q_coeff[5], q_sol[4];
152  FindCoeff(LV_l, LV_l_, LV_b, LV_b_, mt, mt, pxmiss_, pymiss_, q_coeff);
153  int NSol = quartic(q_coeff, q_sol);
154 
155  //loop on all solutions
156  for (int isol = 0; isol < NSol; isol++) {
157  TopRec(LV_l, LV_l_, LV_b, LV_b_, q_sol[isol]);
158  double weight = WeightSolfromShape();
159  if (weight > weightmax) {
160  weightmax =weight;
161  maxLV_n.SetPxPyPzE(LV_n.Px(), LV_n.Py(), LV_n.Pz(), LV_n.E());
162  maxLV_n_.SetPxPyPzE(LV_n_.Px(), LV_n_.Py(), LV_n_.Pz(), LV_n_.E());
163  }
164  }
165  }
167  nuSol.neutrino = reco::LeafCandidate(0, maxLV_n );
168  nuSol.neutrinoBar = reco::LeafCandidate(0, maxLV_n_ );
169  nuSol.weight = weightmax;
170  return nuSol;
171 }
172 
173 void
174 TtFullLepKinSolver::FindCoeff(const TLorentzVector& al,
175  const TLorentzVector& l,
176  const TLorentzVector& b_al,
177  const TLorentzVector& b_l,
178  const double mt,
179  const double mat,
180  const double px_miss,
181  const double py_miss,
182  double* koeficienty)
183 {
184  double E, apom1, apom2, apom3;
185  double k11, k21, k31, k41, cpom1, cpom2, cpom3, l11, l21, l31, l41, l51, l61, k1, k2, k3, k4, k5,k6;
186  double l1, l2, l3, l4, l5, l6, k15, k25, k35, k45;
187 
188  C = -al.Px()-b_al.Px()-l.Px()-b_l.Px() + px_miss;
189  D = -al.Py()-b_al.Py()-l.Py()-b_l.Py() + py_miss;
190 
191  // right side of first two linear equations - missing pT
192 
193  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();
194  F = (sqr(mat)-sqr(mw)-sqr(mb))/(2*b_l.E())-sqr(mw)/(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();
195 
196  m1 = al.Px()/al.E()-b_al.Px()/b_al.E();
197  m2 = al.Py()/al.E()-b_al.Py()/b_al.E();
198  m3 = al.Pz()/al.E()-b_al.Pz()/b_al.E();
199 
200  n1 = l.Px()/l.E()-b_l.Px()/b_l.E();
201  n2 = l.Py()/l.E()-b_l.Py()/b_l.E();
202  n3 = l.Pz()/l.E()-b_l.Pz()/b_l.E();
203 
204  pom = E-m1*C-m2*D;
205  apom1 = sqr(al.Px())-sqr(al.E());
206  apom2 = sqr(al.Py())-sqr(al.E());
207  apom3 = sqr(al.Pz())-sqr(al.E());
208 
209  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);
210  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);
211  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);
212  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));
213  k51 = 1/sqr(al.E())*(apom1*sqr(m3)*sqr(n3)+apom3*sqr(m1)*sqr(n3)-2*al.Px()*al.Pz()*m1*m3*sqr(n3));
214  k61 = 1/sqr(al.E())*(apom2*sqr(m3)*sqr(n3)+apom3*sqr(m2)*sqr(n3)-2*al.Py()*al.Pz()*m2*m3*sqr(n3));
215 
216  cpom1 = sqr(l.Px())-sqr(l.E());
217  cpom2 = sqr(l.Py())-sqr(l.E());
218  cpom3 = sqr(l.Pz())-sqr(l.E());
219 
220  l11 = 1/sqr(l.E())*(pow(mw,4)/4+cpom3*sqr(F)/sqr(n3)+sqr(mw)*l.Pz()*F/n3);
221  l21 = 1/sqr(l.E())*(-2*cpom3*F*m3*n1/n3+sqr(mw)*(l.Px()*m3*n3-l.Pz()*n1*m3)+2*l.Px()*l.Pz()*F*m3);
222  l31 = 1/sqr(l.E())*(-2*cpom3*F*m3*n2/n3+sqr(mw)*(l.Py()*m3*n3-l.Pz()*n2*m3)+2*l.Py()*l.Pz()*F*m3);
223  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));
224  l51 = 1/sqr(l.E())*(cpom1*sqr(m3)*sqr(n3)+cpom3*sqr(n1)*sqr(m3)-2*l.Px()*l.Pz()*n1*n3*sqr(m3));
225  l61 = 1/sqr(l.E())*(cpom2*sqr(m3)*sqr(n3)+cpom3*sqr(n2)*sqr(m3)-2*l.Py()*l.Pz()*n2*n3*sqr(m3));
226 
227  k1 = k11*k61;
228  k2 = k61*k21/k51;
229  k3 = k31;
230  k4 = k41/k51;
231  k5 = k61/k51;
232  k6 = 1;
233 
234  l1 = l11*k61;
235  l2 = l21*k61/k51;
236  l3 = l31;
237  l4 = l41/k51;
238  l5 = l51*k61/(sqr(k51));
239  l6 = l61/k61;
240 
241  k15 = k1*l5-l1*k5;
242  k25 = k2*l5-l2*k5;
243  k35 = k3*l5-l3*k5;
244  k45 = k4*l5-l4*k5;
245 
246  k16 = k1*l6-l1*k6;
247  k26 = k2*l6-l2*k6;
248  k36 = k3*l6-l3*k6;
249  k46 = k4*l6-l4*k6;
250  k56 = k5*l6-l5*k6;
251 
252  koeficienty[0] = k15*sqr(k36)-k35*k36*k16-k56*sqr(k16);
253  koeficienty[1] = 2*k15*k36*k46+k25*sqr(k36)+k35*(-k46*k16-k36*k26)-k45*k36*k16-2*k56*k26*k16;
254  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);
255  koeficienty[3] = k25*sqr(k46)-k35*k46*k56-k45*(k46*k26+k36*k56)-2*sqr(k56)*k26;
256  koeficienty[4] = -k45*k46*k56-pow(k56,3);
257 
258  // normalization of coefficients
259  int moc=(int(log10(fabs(koeficienty[0])))+int(log10(fabs(koeficienty[4]))))/2;
260 
261  koeficienty[0]=koeficienty[0]/TMath::Power(10,moc);
262  koeficienty[1]=koeficienty[1]/TMath::Power(10,moc);
263  koeficienty[2]=koeficienty[2]/TMath::Power(10,moc);
264  koeficienty[3]=koeficienty[3]/TMath::Power(10,moc);
265  koeficienty[4]=koeficienty[4]/TMath::Power(10,moc);
266 }
267 
268 void TtFullLepKinSolver::TopRec(const TLorentzVector& al,
269  const TLorentzVector& l,
270  const TLorentzVector& b_al,
271  const TLorentzVector& b_l,
272  const double sol)
273 {
274  TVector3 t_ttboost;
275  TLorentzVector aux;
276  double pxp, pyp, pzp, pup, pvp, pwp;
277 
278  pxp = sol*(m3*n3/k51);
279  pyp = -(m3*n3/k61)*(k56*pow(sol,2) + k26*sol + k16)/(k36 + k46*sol);
280  pzp = -1/n3*(n1*pxp + n2*pyp - F);
281  pwp = 1/m3*(m1*pxp + m2*pyp + pom);
282  pup = C - pxp;
283  pvp = D - pyp;
284 
285  LV_n_.SetXYZM(pxp, pyp, pzp, 0.0);
286  LV_n.SetXYZM(pup, pvp, pwp, 0.0);
287 
288  LV_t_ = b_l + l + LV_n_;
289  LV_t = b_al + al + LV_n;
290 
291  aux = (LV_t_ + LV_t);
292  t_ttboost = -aux.BoostVector();
293  LV_tt_t_ = LV_t_;
294  LV_tt_t = LV_t;
295  LV_tt_t_.Boost(t_ttboost);
296  LV_tt_t.Boost(t_ttboost);
297 }
298 
299 double
301 {
302  double weight = 1;
303  weight = ((LV_n.E() > genLV_n.E())? genLV_n.E()/LV_n.E(): LV_n.E()/genLV_n.E())
304  *((LV_n_.E() > genLV_n_.E())? genLV_n_.E()/LV_n_.E(): LV_n_.E()/genLV_n_.E());
305  return weight;
306 }
307 
308 double
310 {
311  return EventShape_->Eval(LV_n.E(),LV_n_.E());
312 }
313 
314 int
315 TtFullLepKinSolver::quartic(double *koeficienty, double* koreny) const
316 {
317  double w, b0, b1, b2;
318  double c[4];
319  double d0, d1, h, t, z;
320  double *px;
321 
322  if (koeficienty[4]==0.0)
323  return cubic(koeficienty, koreny);
324  /* quartic problem? */
325  w = koeficienty[3]/(4*koeficienty[4]);
326  /* offset */
327  b2 = -6*sqr(w) + koeficienty[2]/koeficienty[4];
328  /* koeficienty. of shifted polynomial */
329  b1 = (8*sqr(w) - 2*koeficienty[2]/koeficienty[4])*w + koeficienty[1]/koeficienty[4];
330  b0 = ((-3*sqr(w) + koeficienty[2]/koeficienty[4])*w - koeficienty[1]/koeficienty[4])*w + koeficienty[0]/koeficienty[4];
331 
332  c[3] = 1.0;
333  /* cubic resolvent */
334  c[2] = b2;
335  c[1] = -4*b0;
336  c[0] = sqr(b1) - 4*b0*b2;
337 
338  cubic(c, koreny);
339  z = koreny[0];
340  //double z1=1.0,z2=2.0,z3=3.0;
341  //TMath::RootsCubic(c,z1,z2,z3);
342  //if (z2 !=0) z = z2;
343  //if (z1 !=0) z = z1;
344  /* only lowermost root needed */
345 
346  int nreal = 0;
347  px = koreny;
348  t = sqrt(0.25*sqr(z) - b0);
349  for(int i=-1; i<=1; i+=2) {
350  d0 = -0.5*z + i*t;
351  /* coeffs. of quadratic factor */
352  d1 = (t!=0.0)? -i*0.5*b1/t : i*sqrt(-z - b2);
353  h = 0.25*sqr(d1) - d0;
354  if (h>=0.0) {
355  h = sqrt(h);
356  nreal += 2;
357  *px++ = -0.5*d1 - h - w;
358  *px++ = -0.5*d1 + h - w;
359  }
360  }
361 
362  // if (nreal==4) {
363  /* sort results */
364 // if (koreny[2]<koreny[0]) SWAP(koreny[0], koreny[2]);
365 // if (koreny[3]<koreny[1]) SWAP(koreny[1], koreny[3]);
366 // if (koreny[1]<koreny[0]) SWAP(koreny[0], koreny[1]);
367 // if (koreny[3]<koreny[2]) SWAP(koreny[2], koreny[3]);
368 // if (koreny[2]<koreny[1]) SWAP(koreny[1], koreny[2]);
369 // }
370  return nreal;
371 
372 }
373 
374 int
375 TtFullLepKinSolver::cubic(const double *coeffs, double* koreny) const
376 {
377  unsigned nreal;
378  double w, p, q, dis, h, phi;
379 
380  if (coeffs[3]!=0.0) {
381  /* cubic problem? */
382  w = coeffs[2]/(3*coeffs[3]);
383  p = sqr(coeffs[1]/(3*coeffs[3])-sqr(w))*(coeffs[1]/(3*coeffs[3])-sqr(w));
384  q = -0.5*(2*sqr(w)*w-(coeffs[1]*w-coeffs[0])/coeffs[3]);
385  dis = sqr(q)+p;
386  /* discriminant */
387  if (dis<0.0) {
388  /* 3 real solutions */
389  h = q/sqrt(-p);
390  if (h>1.0) h = 1.0;
391  /* confine the argument of */
392  if (h<-1.0) h = -1.0;
393  /* acos to [-1;+1] */
394  phi = acos(h);
395  p = 2*TMath::Power(-p, 1.0/6.0);
396  for(unsigned i=0; i<3; i++)
397  koreny[i] = p*cos((phi+2*i*TMath::Pi())/3.0) - w;
398  if (koreny[1]<koreny[0]) SWAP(koreny[0], koreny[1]);
399  /* sort results */
400  if (koreny[2]<koreny[1]) SWAP(koreny[1], koreny[2]);
401  if (koreny[1]<koreny[0]) SWAP(koreny[0], koreny[1]);
402  nreal = 3;
403  }
404  else {
405  /* only one real solution */
406  dis = sqrt(dis);
407  h = TMath::Power(fabs(q+dis), 1.0/3.0);
408  p = TMath::Power(fabs(q-dis), 1.0/3.0);
409  koreny[0] = ((q+dis>0.0)? h : -h) + ((q-dis>0.0)? p : -p) - w;
410  nreal = 1;
411  }
412 
413  /* Perform one step of a Newton iteration in order to minimize
414  round-off errors */
415  for(unsigned i=0; i<nreal; i++) {
416  h = coeffs[1] + koreny[i] * (2 * coeffs[2] + 3 * koreny[i] * coeffs[3]);
417  if (h != 0.0)
418  koreny[i] -= (coeffs[0] + koreny[i] * (coeffs[1] + koreny[i] * (coeffs[2] + koreny[i] * coeffs[3])))/h;
419  }
420  }
421 
422  else if (coeffs[2]!=0.0) {
423  /* quadratic problem? */
424  p = 0.5*coeffs[1]/coeffs[2];
425  dis = sqr(p) - coeffs[0]/coeffs[2];
426  if (dis>=0.0) {
427  /* two real solutions */
428  dis = sqrt(dis);
429  koreny[0] = -p - dis;
430  koreny[1] = -p + dis;
431  nreal = 2;
432  }
433  else
434  /* no real solution */
435  nreal = 0;
436  }
437 
438  else if (coeffs[1]!=0.0) {
439  /* linear problem? */
440  koreny[0] = -coeffs[0]/coeffs[1];
441  nreal = 1;
442  }
443 
444  else
445  /* no equation */
446  nreal = 0;
447 
448  return nreal;
449 }
450 
451 
452 void
453 TtFullLepKinSolver::SWAP(double& realone, double& realtwo) const
454 {
455  if (realtwo < realone) {
456  double aux = realtwo;
457  realtwo = realone;
458  realone = aux;
459  }
460 }
const double Pi
pat::Jet getCalJetBbar() const
void setRecWeightMax(double wgt)
tuple t
Definition: tree.py:139
int i
Definition: DBlmapReader.cc:9
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_)
const double w
Definition: UKUtility.cc:23
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()
const reco::GenParticle * getGenN() const
pat::Electron getElectronm() const
std::string getWmDecay() const
std::string getWpDecay() const
TLorentzVector genLV_n
provisional
const double topmass_step
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
void SWAP(double &realone, double &realtwo) const
virtual double energy() const
energy
const reco::GenParticle * getGenNbar() const
TtDilepEvtSolution addKinSolInfo(TtDilepEvtSolution *asol)
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
TF2 * EventShape_
Event shape.
double sqr(const double x) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
TLorentzVector genLV_n_
pat::Electron getElectronp() const
void TopRec(const TLorentzVector &al, const TLorentzVector &l, const TLorentzVector &b_al, const TLorentzVector &b_l, const double sol)
pat::Tau getTaup() const
~TtFullLepKinSolver()
destructor
pat::Muon getMuonm() const
double WeightSolfromShape() const
use the parametrized event shape to obtain the solution weight.
virtual double px() const
x coordinate of momentum vector
pat::Tau getTaum() const
pat::Jet getCalJetB() const
double b
Definition: hdecay.h:120
virtual double pz() const
z coordinate of momentum vector
TLorentzVector LV_tt_t
int weight
Definition: histoStyle.py:50
pat::Muon getMuonp() const
const double topmass_begin
TLorentzVector LV_tt_t_
void setRecTopMass(double mass)
virtual double py() const
y coordinate of momentum vector
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40