Belle II Software development
EvtBSemiTauonic2HDMType2.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#include <string>
14
15#include "framework/logging/Logger.h"
16
17#include "generators/evtgen/EvtGenModelRegister.h"
18#include "generators/evtgen/models/EvtBSemiTauonic2HDMType2.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 EvtBSemiTauonic2HDMType2::EvtBSemiTauonic2HDMType2(): 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_2HDMTYPE2";
43 }
44
45
46
48 {
49 return new EvtBSemiTauonic2HDMType2;
50 }
51
52
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("EvtBSemiTauonic2HDMType2::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 const double m_tau = EvtPDL::getMeanMass(getDaug(1));
90
91 EvtSpinType::spintype d1type = EvtPDL::getSpinType(getDaug(0));
92
93 if (d1type == EvtSpinType::VECTOR) {
94 B2INFO("EvtBSemiTauonic2HDMType2::init()>> Initializing for decays to a vector type meson ");
95 checkNArg(7, 8);
96
97 const double rhoa12 = getArg(0);
98 const double R11 = getArg(1);
99 const double R21 = getArg(2);
100 const double aR3 = getArg(3);
101
102 const double rho12 = 0; // not used for D*taunu decays
103 const double aS1 = 0; // not used for D*taunu decays
104
105 EvtComplex Coeffs[5];
106
107 const double m_b = getArg(4);
108 const double m_c = getArg(5);
109 const double tanBetaOverMH = getArg(6);
110 B2INFO("tan(beta)/m_H+ = " << tanBetaOverMH);
111
112 Coeffs[0] = 0; // CV1
113 Coeffs[1] = 0; // CV2
114 Coeffs[2] = -m_b * m_tau * tanBetaOverMH * tanBetaOverMH; // CS1
115 Coeffs[3] = 0; // CS2, neglected by default
116 Coeffs[4] = 0; // CT
117
118 if (getNArg() == 8) {
119 // if m_H is explicitly given
120 const double m_H = getArg(7);
121 Coeffs[3] = -m_c * m_tau / m_H / m_H; // CS2
122 }
123
124 m_CalcHelAmp = new EvtBSemiTauonicHelicityAmplitudeCalculator(rho12, rhoa12, R11, R21, aS1, aR3, m_b, m_c,
125 Coeffs[0], Coeffs[1], Coeffs[2], Coeffs[3], Coeffs[4],
126 EvtPDL::getMeanMass(getParentId()),
127 -1, /*dummy for Dmass*/
128 EvtPDL::getMeanMass(getDaug(0)));
130 } else if (d1type == EvtSpinType::SCALAR) {
131 B2INFO("EvtBSemiTauonic2HDMType2::init()>> Initializing for decays to a scalar type meson ");
132 checkNArg(5, 6);
133
134 const double rho12 = getArg(0);
135 const double aS1 = getArg(1);
136
137 const double rhoa12 = 0; // not used for Dtaunu decays
138 const double R11 = 0; // not used for Dtaunu decays
139 const double R21 = 0; // not used for Dtaunu decays
140 const double aR3 = 0; // not used for Dtaunu decays
141
142 EvtComplex Coeffs[5];
143
144 const double m_b = getArg(2);
145 const double m_c = getArg(3);
146 const double tanBetaOverMH = getArg(4);
147 B2INFO("tan(beta)/m_H+ = " << tanBetaOverMH);
148
149 Coeffs[0] = 0; // CV1
150 Coeffs[1] = 0; // CV2
151 Coeffs[2] = -m_b * m_tau * tanBetaOverMH * tanBetaOverMH; // CS1
152 Coeffs[3] = 0; // CS2, neglected by default
153 Coeffs[4] = 0; // CT
154
155 if (getNArg() == 6) {
156 // if m_H is explicitly given
157 const double m_H = getArg(5);
158 Coeffs[3] = -m_c * m_tau / m_H / m_H; // CS2
159 }
160
161 m_CalcHelAmp = new EvtBSemiTauonicHelicityAmplitudeCalculator(rho12, rhoa12, R11, R21, aS1, aR3, m_b, m_c,
162 Coeffs[0], Coeffs[1], Coeffs[2], Coeffs[3], Coeffs[4],
163 EvtPDL::getMeanMass(getParentId()),
164 EvtPDL::getMeanMass(getDaug(0)),
165 -1 /*dummy for D*mass*/);
167 } else {
168 B2ERROR("BSemiTauonic2HDMType2 model handles only scalar and vector meson daughters.");
169// EvtGenReport(EVTGEN_ERROR, "EvtGen") << "BSemiTauonic2HDMType2 model handles only scalar and vector meson daughters."<<std::endl;
170 ::abort();
171 }
172 }
174} // Belle 2 Namespace
The EvtGen model of semi-tauonic B decays including the charged higgs effect of 2HDM Type2 based on [...
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.
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...
void init()
The function initializes the decay.
EvtBSemiTauonic2HDMType2()
The default constructor
virtual ~EvtBSemiTauonic2HDMType2()
The destructor
EvtDecayBase * clone()
The function makes a copy of an EvtBSemiTauonic2HDMType2 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.
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.