Belle II Software  release-06-02-00
SpecificKinematicVariables.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 // Own include
10 #include <analysis/variables/SpecificKinematicVariables.h>
11 #include <analysis/utility/PCmsLabTransform.h>
12 #include <analysis/variables/Variables.h>
13 
14 // dataobjects
15 #include <analysis/dataobjects/Particle.h>
16 
17 // framework aux
18 #include <framework/logging/Logger.h>
19 
20 #include <TRandom.h>
21 #include <TMath.h>
22 #include <iostream>
23 using namespace std;
24 
25 namespace Belle2 {
30  namespace Variable {
31 
32 
33  double REC_q2BhSimple(const Particle* particle)
34  {
35  // calculates q^2 = (p_B - p_h) in decays of B -> h_1 .. h_n ell nu_ell,
36  // where p_h = Sum_i^n p_h_i is the 4-momentum of hadrons in the final
37  // state. The calculation is performed in the CMS system, where B-meson
38  // is assumed to be at rest p_B = (m_B, 0).
39 
40  TLorentzVector hadron4vec;
41 
42  unsigned n = particle->getNDaughters();
43 
44  if (n < 1)
45  return std::numeric_limits<float>::quiet_NaN();
46 
47  for (unsigned i = 0; i < n; i++) {
48  int absPDG = abs(particle->getDaughter(i)->getPDGCode());
49  if (absPDG == Const::electron.getPDGCode() || absPDG == Const::muon.getPDGCode() || absPDG == 15)
50  continue;
51 
52  hadron4vec += particle->getDaughter(i)->get4Vector();
53  }
54 
55  // boost to CMS
56  PCmsLabTransform T;
57  TLorentzVector phCMS = T.rotateLabToCms() * hadron4vec;
58  TLorentzVector pBCMS;
59  pBCMS.SetXYZM(0.0, 0.0, 0.0, particle->getPDGMass());
60  return (pBCMS - phCMS).Mag2();
61  }
62 
63  double REC_q2Bh(const Particle* particle)
64  {
65  // calculates q^2 = (p_B - p_h) in decays of B -> h_1 .. h_n ell nu_ell,
66  // where p_h = Sum_i^n p_h_i is the 4-momentum of hadrons in the final
67  // state. The calculation is performed in the CMS system,
68  // with a weighter average in a cone around the true B direction
69 
70  TLorentzVector hadron4vec;
71 
72  unsigned n = particle->getNDaughters();
73 
74  if (n < 1)
75  return std::numeric_limits<float>::quiet_NaN();
76 
77  for (unsigned i = 0; i < n; i++) {
78  int absPDG = abs(particle->getDaughter(i)->getPDGCode());
79  if (absPDG == Const::electron.getPDGCode() || absPDG == Const::muon.getPDGCode() || absPDG == 15)
80  continue;
81 
82  hadron4vec += particle->getDaughter(i)->get4Vector();
83  }
84 
85  // boost to CMS
86  PCmsLabTransform T;
87  TLorentzVector had_cm = T.rotateLabToCms() * hadron4vec;
88  TLorentzVector Y_cm = T.rotateLabToCms() * particle->get4Vector();
89 
90 
91  double bmass = particle->getPDGMass();
92  // B theta angle
93  double cos_cone_angle = Variable::cosThetaBetweenParticleAndNominalB(particle);
94  if (abs(cos_cone_angle) > 1) {
95  //makes no sense in this case, return simple value
96  double q2 = Variable::REC_q2BhSimple(particle);
97  if (q2 < 0) {
98  return 0;
99  }
100  return q2;
101  }
102  double thetaBY = TMath::ACos(cos_cone_angle);
103  const double E_B = T.getCMSEnergy() / 2.0;
104  const double p_B = sqrt(E_B * E_B - bmass * bmass);
105 
106  double phi_start = gRandom->Uniform(0, TMath::Pi() / 2);
107 
108  double q2 = 0;
109  double denom = 0;
110 
111  TVector3 zHatY(Y_cm.Vect().Unit());
112  TVector3 yHatY((had_cm.Vect().Cross(zHatY)).Unit());
113  TVector3 xHatY(yHatY.Cross(zHatY));
114 
115  double phi = phi_start;
116  double dphi = TMath::Pi() / 2;
117  TLorentzVector static_B(0, 0, 0, bmass);
118 
119  for (int around_the_cone = 0; around_the_cone < 4; around_the_cone++) {
120  //Define the momentum of B in the Y rest frame using the angles thetaBY and the
121  //current phi in the loop.
122 
123  //The purpose of these lines is to calculate the B momentum in the BY cone in the coordinate system previously developed.
124  double B0_px_Y_frame = p_B * TMath::Sin(thetaBY) * TMath::Cos(phi);
125  double B0_py_Y_frame = p_B * TMath::Sin(thetaBY) * TMath::Sin(phi);
126  double B0_pz_Y_frame = p_B * TMath::Cos(thetaBY);
127 
128  //The 3 components of the guess for the B0 p direction at the current phi are
129  //calculated by scaling the 3 basis vectors computed before by the corresponding B
130  //momentum in that direction.
131 
132  TVector3 B0_px_Dframe(xHatY);
133  B0_px_Dframe *= B0_px_Y_frame;
134  TVector3 B0_py_Dframe(yHatY);
135  B0_py_Dframe *= B0_py_Y_frame;
136  TVector3 B0_pz_Dframe(zHatY);
137  B0_pz_Dframe *= B0_pz_Y_frame;
138  //Construct the B0 p 3-vector with the current phi by summing the 3 components.
139 
140  TVector3 B0_p3_Dframe(B0_px_Dframe + B0_py_Dframe + B0_pz_Dframe);
141  TLorentzVector B0_p4_Dframe(B0_p3_Dframe, E_B);
142 
143  //This is the polar angle of B0.
144 
145  double cosThetaB = B0_p3_Dframe.CosTheta();
146  double sinThetaB2 = (1 - cosThetaB * cosThetaB);
147 
148  //The weight is given by the sin squared of such angle.
149  double wt = sinThetaB2;
150 
151  //Boost the hadronic daughter to the computed B0 rest frame.
152  // In that frame, q2 can simply be calculated by subtracting the hadron 4-momentum from the momentum of a static B.
153  TLorentzVector had_B0(had_cm);
154  had_B0.Boost(-B0_p4_Dframe.BoostVector());
155  q2 += wt * ((static_B - had_B0).M2());
156  denom += wt;
157  phi += dphi;
158  }
159 
160  q2 /= denom;
161  if (q2 < 0) {
162  q2 = 0.0;
163  }
164  return q2;
165  }
166 
167  double REC_MissM2(const Particle* particle)
168  {
169  PCmsLabTransform T;
170  TLorentzVector rec4vecLAB = particle->get4Vector();
171  TLorentzVector rec4vec = T.rotateLabToCms() * rec4vecLAB;
172 
173  TLorentzVector miss4vec;
174  double E_beam_cms = T.getCMSEnergy() / 2.0;
175 
176  miss4vec.SetVect(-rec4vec.Vect());
177  miss4vec.SetE(E_beam_cms - rec4vec.Energy());
178 
179  return miss4vec.Mag2();
180  }
181 
182 
183  VARIABLE_GROUP("Specific kinematic variables");
184 
185  REGISTER_VARIABLE("recQ2BhSimple", REC_q2BhSimple,
186  "Returns the momentum transfer squared, :math:`q^2`, calculated in CMS as :math:`q^2 = (p_B - p_h)^2`, \n"
187  "where p_h is the CMS momentum of all hadrons in the decay :math:`B \\to H_1 ... H_n \\ell \\nu_\\ell`.\n"
188  "The B meson momentum in CMS is assumed to be 0.");
189 
190  REGISTER_VARIABLE("recQ2Bh", REC_q2Bh,
191  "Returns the momentum transfer squared, :math:`q^2`, calculated in CMS as :math:`q^2 = (p_B - p_h)^2`, \n"
192  "where p_h is the CMS momentum of all hadrons in the decay :math:`B \\to H_1\\dots H_n \\ell \\nu_\\ell`.\n"
193  "This calculation uses a weighted average of the B meson around the reco B cone. \n"
194  "Based on diamond frame calculation of :math:`q^2` following the idea presented in https://www.osti.gov/biblio/1442697 \n"
195  "It will switch to use of :b2:var:`recQ2BhSimple` if absolute of :b2:var:`cosThetaBetweenParticleAndNominalB` > 1.");
196 
197  REGISTER_VARIABLE("recMissM2", REC_MissM2,
198  "Returns the invariant mass squared of the missing momentum calculated assumings the"
199  "reco B is at rest and calculating the neutrino (missing) momentum from :math:`p_\\nu = p_B - p_\\mathrm{had} - p_\\mathrm{lep}`");
200 
201 
202 
203  }
205 }
Abstract base class for different kinds of events.