10#include <analysis/variables/SpecificKinematicVariables.h>
13#include <analysis/VariableManager/Manager.h>
15#include <analysis/utility/PCmsLabTransform.h>
16#include <analysis/variables/Variables.h>
19#include <analysis/dataobjects/Particle.h>
23#include <Math/Vector3D.h>
24#include <Math/Vector4D.h>
26using namespace ROOT::Math;
37 double REC_q2BhSimple(
const Particle* particle)
44 PxPyPzEVector hadron4vec;
46 unsigned n = particle->getNDaughters();
51 for (
unsigned i = 0; i < n; i++) {
52 int absPDG = abs(particle->getDaughter(i)->getPDGCode());
56 hadron4vec += particle->getDaughter(i)->get4Vector();
61 PxPyPzEVector phCMS = T.rotateLabToCms() * hadron4vec;
62 PxPyPzMVector pBCMS(0.0, 0.0, 0.0, particle->getPDGMass());
63 return (pBCMS - phCMS).mag2();
66 double REC_q2Bh(
const Particle* particle)
73 PxPyPzEVector hadron4vec;
75 unsigned n = particle->getNDaughters();
80 for (
unsigned i = 0; i < n; i++) {
81 int absPDG = abs(particle->getDaughter(i)->getPDGCode());
85 hadron4vec += particle->getDaughter(i)->get4Vector();
90 PxPyPzEVector had_cm = T.rotateLabToCms() * hadron4vec;
91 PxPyPzEVector Y_cm = T.rotateLabToCms() * particle->get4Vector();
94 double bmass = particle->getPDGMass();
96 double cos_cone_angle = Variable::cosThetaBetweenParticleAndNominalB(particle);
97 if (abs(cos_cone_angle) > 1) {
99 double q2 = Variable::REC_q2BhSimple(particle);
105 double thetaBY = TMath::ACos(cos_cone_angle);
106 const double E_B = T.getCMSEnergy() / 2.0;
107 const double p_B =
sqrt(E_B * E_B - bmass * bmass);
109 double phi_start = gRandom->Uniform(0, TMath::Pi() / 2);
114 XYZVector zHatY(Y_cm.Vect().Unit());
115 XYZVector yHatY((had_cm.Vect().Cross(zHatY)).Unit());
116 XYZVector xHatY(yHatY.Cross(zHatY));
118 double phi = phi_start;
119 double dphi = TMath::Pi() / 2;
120 PxPyPzEVector static_B(0, 0, 0, bmass);
122 for (
int around_the_cone = 0; around_the_cone < 4; around_the_cone++) {
127 double B0_px_Y_frame = p_B * TMath::Sin(thetaBY) * TMath::Cos(phi);
128 double B0_py_Y_frame = p_B * TMath::Sin(thetaBY) * TMath::Sin(phi);
129 double B0_pz_Y_frame = p_B * TMath::Cos(thetaBY);
135 XYZVector B0_px_Dframe(xHatY);
136 B0_px_Dframe *= B0_px_Y_frame;
137 XYZVector B0_py_Dframe(yHatY);
138 B0_py_Dframe *= B0_py_Y_frame;
139 XYZVector B0_pz_Dframe(zHatY);
140 B0_pz_Dframe *= B0_pz_Y_frame;
143 XYZVector B0_p3_Dframe(B0_px_Dframe + B0_py_Dframe + B0_pz_Dframe);
144 PxPyPzEVector B0_p4_Dframe(B0_p3_Dframe.X(), B0_p3_Dframe.Y(), B0_p3_Dframe.Z(), E_B);
148 double cosThetaB = cos(B0_p3_Dframe.Theta());
149 double sinThetaB2 = (1 - cosThetaB * cosThetaB);
152 double wt = sinThetaB2;
156 PxPyPzEVector had_B0(had_cm);
157 had_B0 = Boost(B0_p4_Dframe.BoostToCM()) * had_B0;
158 q2 += wt * ((static_B - had_B0).M2());
170 double REC_MissM2(
const Particle* particle)
173 PxPyPzEVector rec4vecLAB = particle->get4Vector();
174 PxPyPzEVector rec4vec = T.rotateLabToCms() * rec4vecLAB;
176 double E_beam_cms = T.getCMSEnergy() / 2.0;
178 PxPyPzEVector miss4vec(-rec4vec.px(), -rec4vec.py(), -rec4vec.pz(), E_beam_cms - rec4vec.E());
180 return miss4vec.mag2();
184 VARIABLE_GROUP(
"Specific kinematic variables");
186 REGISTER_VARIABLE(
"recQ2BhSimple", REC_q2BhSimple,
187 "Returns the momentum transfer squared, :math:`q^2`, calculated in CMS as :math:`q^2 = (p_B - p_h)^2`, \n"
188 "where p_h is the CMS momentum of all hadrons in the decay :math:`B \\to H_1 ... H_n \\ell \\nu_\\ell`.\n"
189 "The B meson momentum in CMS is assumed to be 0.\n\n",
":math:`[\\text{GeV}/\\text{c}]^2`");
191 REGISTER_VARIABLE(
"recQ2Bh", REC_q2Bh,
192 "Returns the momentum transfer squared, :math:`q^2`, calculated in CMS as :math:`q^2 = (p_B - p_h)^2`, \n"
193 "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"
194 "This calculation uses a weighted average of the B meson around the reco B cone. \n"
195 "Based on diamond frame calculation of :math:`q^2` following the idea presented in https://www.osti.gov/biblio/1442697 \n"
196 "It will switch to use of :b2:var:`recQ2BhSimple` if absolute of :b2:var:`cosThetaBetweenParticleAndNominalB` > 1.\n\n",
197 ":math:`[\\text{GeV}/\\text{c}]^2`");
199 REGISTER_VARIABLE(
"recMissM2", REC_MissM2, R
"DOC(
200 Returns the invariant mass squared of the missing momentum calculated assumings the
201 reco B is at rest and calculating the neutrino (missing) momentum from :math:`p_\nu = p_B - p_{\rm had} - p_{\rm lep}`
203 )DOC", ":math:`[\\text{GeV}/\\text{c}^2]^2`");
static const ChargedStable muon
muon particle
static const double doubleNaN
quiet_NaN
static const ChargedStable electron
electron particle
double sqrt(double a)
sqrt for double
Abstract base class for different kinds of events.