Belle II Software development
EvtBSemiTauonic.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#include <stdlib.h>
10#include "EvtGenBase/EvtParticle.hh"
11#include "EvtGenBase/EvtPDL.hh"
12
13
14#include <string>
15
16#include "framework/logging/Logger.h"
17#include <generators/evtgen/EvtGenModelRegister.h>
18#include "generators/evtgen/models/EvtBSemiTauonic.h"
19#include "generators/evtgen/models/EvtBSemiTauonicHelicityAmplitudeCalculator.h"
20#include "generators/evtgen/models/EvtBSemiTauonicVectorMesonAmplitude.h"
21#include "generators/evtgen/models/EvtBSemiTauonicScalarMesonAmplitude.h"
22
23namespace Belle2 {
31
32 EvtBSemiTauonic::EvtBSemiTauonic(): m_CalcHelAmp(0), m_CalcAmp(0) {}
33
35 {
36 if (m_CalcHelAmp)delete m_CalcHelAmp;
37 if (m_CalcAmp)delete m_CalcAmp;
38 }
39
41 {
42 return "BSTD";
43 }
44
45
46
47 EvtDecayBase* EvtBSemiTauonic::clone()
48 {
49 return new EvtBSemiTauonic;
50 }
51
52
53 void EvtBSemiTauonic::decay(EvtParticle* p)
54 {
55 p->initializePhaseSpace(getNDaug(), getDaugs());
56 m_CalcAmp->CalcAmp(p, _amp2, m_CalcHelAmp);
57 }
58
60 {
61 EvtId parnum, mesnum, lnum, nunum;
62
63 parnum = getParentId();
64 mesnum = getDaug(0);
65 lnum = getDaug(1);
66 nunum = getDaug(2);
67
68 double maxprob = m_CalcAmp->CalcMaxProb(parnum, mesnum,
69 lnum, nunum,
71
72 B2INFO("EvtBSemiTauonic::initProbMax()>> maxprob: " << maxprob);
73
74 setProbMax(maxprob);
75 }
76
77
79 {
80 checkNDaug(3);
81
82 //We expect the parent to be a scalar
83 //and the daughters to be X lepton neutrino
84 checkSpinParent(EvtSpinType::SCALAR);
85
86 checkSpinDaughter(1, EvtSpinType::DIRAC);
87 checkSpinDaughter(2, EvtSpinType::NEUTRINO);
88
89 EvtSpinType::spintype d1type = EvtPDL::getSpinType(getDaug(0));
90
91 if (d1type == EvtSpinType::VECTOR) {
92 B2INFO("EvtBSemiTauonic::init()>> Initializing for decays to a vector type meson ");
93 checkNArg(16);
94 const double rhoa12 = getArg(0);
95 const double R11 = getArg(1);
96 const double R21 = getArg(2);
97 const double aR3 = getArg(3);
98
99 const double rho12 = 0; // not used for D*taunu decays
100 const double aS1 = 0; // not used for D*taunu decays
101
102 const double m_b = getArg(4);
103 const double m_c = getArg(5);
104
105 EvtComplex Coeffs[5];
106 Coeffs[0] = EvtComplex(getArg(6) * cos(getArg(7)), getArg(6) * sin(getArg(7)));
107 Coeffs[1] = EvtComplex(getArg(8) * cos(getArg(9)), getArg(8) * sin(getArg(9)));
108 Coeffs[2] = EvtComplex(getArg(10) * cos(getArg(11)), getArg(10) * sin(getArg(11)));
109 Coeffs[3] = EvtComplex(getArg(12) * cos(getArg(13)), getArg(12) * sin(getArg(13)));
110 Coeffs[4] = EvtComplex(getArg(14) * cos(getArg(15)), getArg(14) * sin(getArg(15)));
111
112 m_CalcHelAmp = new EvtBSemiTauonicHelicityAmplitudeCalculator(rho12, rhoa12, R11, R21, aS1, aR3, m_b, m_c,
113 Coeffs[0], Coeffs[1], Coeffs[2], Coeffs[3], Coeffs[4],
114 EvtPDL::getMeanMass(getParentId()),
115 -1, /*dummy for Dmass*/
116 EvtPDL::getMeanMass(getDaug(0)));
118 } else if (d1type == EvtSpinType::SCALAR) {
119 B2INFO("EvtBSemiTauonic::init()>> Initializing for decays to a scalar type meson ");
120 checkNArg(14);
121 const double rho12 = getArg(0);
122 const double aS1 = getArg(1);
123
124 const double rhoa12 = 0; // not used for Dtaunu decays
125 const double R11 = 0; // not used for Dtaunu decays
126 const double R21 = 0; // not used for Dtaunu decays
127 const double aR3 = 0; // not used for Dtaunu decays
128
129 const double m_b = getArg(2);
130 const double m_c = getArg(3);
131
132 EvtComplex Coeffs[5];
133 Coeffs[0] = EvtComplex(getArg(4) * cos(getArg(5)), getArg(4) * sin(getArg(5)));
134 Coeffs[1] = EvtComplex(getArg(6) * cos(getArg(7)), getArg(6) * sin(getArg(7)));
135 Coeffs[2] = EvtComplex(getArg(8) * cos(getArg(9)), getArg(8) * sin(getArg(9)));
136 Coeffs[3] = EvtComplex(getArg(10) * cos(getArg(11)), getArg(10) * sin(getArg(11)));
137 Coeffs[4] = EvtComplex(getArg(12) * cos(getArg(13)), getArg(12) * sin(getArg(13)));
138
139 m_CalcHelAmp = new EvtBSemiTauonicHelicityAmplitudeCalculator(rho12, rhoa12, R11, R21, aS1, aR3, m_b, m_c,
140 Coeffs[0], Coeffs[1], Coeffs[2], Coeffs[3], Coeffs[4],
141 EvtPDL::getMeanMass(getParentId()),
142 EvtPDL::getMeanMass(getDaug(0)),
143 -1 /*dummy for D*mass*/);
145 } else {
146 B2ERROR("BSemiTauonic model handles only scalar and vector meson daughters.");
147// EvtGenReport(EVTGEN_ERROR, "EvtGen") << "BSemiTauonic model handles only scalar and vector meson daughters." << std::endl;
148 ::abort();
149 }
150 }
152} // Belle 2 Namespace
virtual void CalcAmp(EvtParticle *parent, EvtAmp &amp, EvtBSemiTauonicHelicityAmplitudeCalculator *HelicityAmplitudeCalculator)=0
The function calculates the spin dependent amplitude.
The class calculates the helicity amplitude of semi-tauonic B decays including new physics effects ba...
The class calculates the spin dependent amplitudes of B->Dtaunu decays for the BSemiTauonic model bas...
The class calculates the spin dependent amplitudes of B->D*taunu decays for the BSemiTauonic model ba...
The EvtGen model of semi-tauonic B decays including new physics effects based on [M.
EvtBSemiTauonicAmplitude * m_CalcAmp
A pointer to the spin-dependent amplitude calculator specific to the spin type of the daughter meson.
EvtBSemiTauonicHelicityAmplitudeCalculator * m_CalcHelAmp
A pointer to the helicity amplitude calculator.
void init()
The function initializes the decay.
EvtDecayBase * clone()
The function makes a copy of an EvtBSTD object.
void initProbMax()
The function sets the maximum value of the probability.
std::string getName()
The function returns the model name.
#define B2_EVTGEN_REGISTER_MODEL(classname)
Class to register B2_EVTGEN_REGISTER_MODEL.
EvtBSemiTauonic()
The default constructor
virtual ~EvtBSemiTauonic()
The destructor
void decay(EvtParticle *p)
The function evaluates the decay amplitude of the parent particle.
double CalcMaxProb(EvtId parent, EvtId meson, EvtId lepton, EvtId nudaug, EvtBSemiTauonicHelicityAmplitudeCalculator *HelicityAmplitudeCalculator)
The function calculates the maximum probability.
Abstract base class for different kinds of events.