Belle II Software  release-05-01-25
EvtBSemiTauonicHelicityAmplitudeCalculator.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Koji Hara *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #pragma once
11 
12 #include "EvtGenBase/EvtComplex.hh"
13 
14 namespace Belle2 {
30  class EvtBSemiTauonicHelicityAmplitudeCalculator {
31 
32  public:
33 
38 
59  EvtBSemiTauonicHelicityAmplitudeCalculator(const double rho12, const double rhoA12, const double ffR11, const double ffR21,
60  const double AS1, const double AR3,
61  const double bottomMass, const double charmMass,
62  const EvtComplex& CV1, const EvtComplex& CV2, const EvtComplex& CS1, const EvtComplex& CS2, const EvtComplex& CT,
63  const double parentMass, const double DMass, const double DstarMass);
64 
67 
76  EvtComplex helAmp(double mtau, int tauhel, int Dhel, double w, double costau) const;
77 
78  public:
79 
93  EvtComplex helAmp(const EvtComplex& CV1, const EvtComplex& CV2, const EvtComplex& CS1, const EvtComplex& CS2, const EvtComplex& CT,
94  double mtau, int tauhel, int Dhel, double w, double costau) const;
95 
104  double Lep(const double mtau, int tauhel, int whel, double q2, double costau) const; // Vector
105 
113  double Lep(const double mtau, int tauhel, double q2, double costau) const; // Scalar
114 
125  double Lep(const double mtau, int tauhel, int whel1, int whel2, double q2, double costau) const; // Tensor
126 
133  double HadV1(int Dhel, int whel, double w) const; // V-A
134 
141  double HadV2(int Dhel, int whel, double w) const; // V+A
142 
149  double HadS1(int Dhel, double w) const; // S+P
150 
157  double HadS2(int Dhel, double w) const; // S-P
158 
167  double HadT(int Dhel, int whel1, int whel2, double w) const; // Tensor
168 
178  double helampSM(double mtau, int tauhel, int Dhel, double w, double costau) const; // SM
179 
190  double helampV1(double mtau, int tauhel, int Dhel, double w, double costau) const; // V-A
191 
202  double helampV2(double mtau, int tauhel, int Dhel, double w, double costau) const; // V+A
203 
214  double helampS1(double mtau, int tauhel, int Dhel, double w, double costau) const; // S+P
215 
226  double helampS2(double mtau, int tauhel, int Dhel, double w, double costau) const; // S-P
227 
238  double helampT(double mtau, int tauhel, int Dhel, double w, double costau) const; // Tensor
239 
244  double eta(int whel) const {return (whel == 2) ? -1 : 1;}
245 
246 
248 
249  // Vector and axial-vector form factors in terms of CLN form factors
250 
252  double hp(double w) const;
253 
255  double hm(double w) const;
256 
258  double hA1(double w) const;
259 
261  double hV(double w) const;
262 
264  double hA2(double w) const;
265 
267  double hA3(double w) const;
268 
269  // Scalar form factors in terms of V and A form factors
270 
272  double hS(double w) const;
273 
275  double hP(double w) const;
276 
277  // Tensor form factors in terms of V and A form factors
278 
280  double hT(double w) const;
281 
283  double hT1(double w) const;
284 
286  double hT2(double w) const;
287 
289  double hT3(double w) const;
290 
295  double z(double w) const;
296 
301  double ffV1(double w) const;
302 
307  double ffS1(double w) const;
308 
313  double ffA1(double w) const;
314 
319  double ffR1(double w) const;
320 
325  double ffR2(double w) const;
326 
331  double ffR3(double w) const;
332 
337  double ffV11() const {return 1.;} // cancels in R(D)
338 
343  double ffA11() const {return 1.;} // cancels in R(Ds)
344 
349  double dS1(double w) const;
350 
355  double dR3(double w) const;
356 
360  double aS1() const {return getAS1();} // {return 1.;}
361 
365  double aR3() const {return getAR3();} // {return 1.;}
366 
371  double mD(int Dhel) const;
372 
377  double r(int Dhel) const {return mD(Dhel) / m_mB;} // meson mass ratio
378 
379 
382  double rq() const {return m_mCharm / m_mBottom;} // quark mass ratio
383 
389  double v(double mtau, double q2) const;
390 
396  double q2(int Dhel, double w) const;
397 
403  double qh2(int Dhel, double w) const;
404 
405  // range of q2 and w
406 
411  double q2min(double mtau) const {return mtau * mtau;}
412 
417  double q2max(int Dhel) const {return (m_mB - mD(Dhel)) * (m_mB - mD(Dhel));}
418 
424  double wfunc(int Dhel, double q2) const {return (1. + r(Dhel) * r(Dhel) - q2 / m_mB / m_mB) / 2. / r(Dhel);}
425 
429  double wmin() const {return 1.;};
430 
436  double wmax(double mtau, int Dhel) const {return wfunc(Dhel, q2min(mtau));}
437 
441  double getRho12() const {return m_rho12;}
442 
444  double getRhoA12() const {return m_rhoA12;}
445 
447  double getR11() const {return m_ffR11;}
448 
450  double getR21() const {return m_ffR21;}
451 
453  double getAS1() const {return m_aS1;}
454 
456  double getAR3() const {return m_aR3;}
457 
459  double getMB() const {return m_mB;}
460 
462  double getMD() const {return m_mD;}
463 
465  double getMDst() const {return m_mDst;}
466 
468  double getMBottom() const {return m_mBottom;}
469 
471  double getMCharm() const {return m_mCharm;}
472 
474  EvtComplex getCV1() const {return m_CV1;}
475 
477  EvtComplex getCV2() const {return m_CV2;}
478 
480  EvtComplex getCS1() const {return m_CS1;}
481 
483  EvtComplex getCS2() const {return m_CS2;}
484 
486  EvtComplex getCT() const {return m_CT;}
487 
489  void setRho12(double v) {m_rho12 = v;}
490 
492  void setRhoA12(double v) {m_rhoA12 = v;}
493 
495  void setR11(double v) {m_ffR11 = v;}
496 
498  void setR21(double v) {m_ffR21 = v;}
499 
501  void setAS1(double v) {m_aS1 = v;}
502 
504  void setAR3(double v) {m_aR3 = v;}
505 
507  void setMB(double m) {m_mB = m;}
508 
510  void setMD(double m) {m_mD = m;}
511 
513  void setMDst(double m) {m_mDst = m;}
514 
516  void setMBottom(double m) {m_mBottom = m;}
517 
519  void setMCharm(double m) {m_mCharm = m;}
520 
522  void setCV1(const EvtComplex& v) {m_CV1 = v;}
523 
525  void setCV2(const EvtComplex& v) {m_CV2 = v;}
526 
528  void setCS1(const EvtComplex& v) {m_CS1 = v;}
529 
531  void setCS2(const EvtComplex& v) {m_CS2 = v;}
532 
534  void setCT(const EvtComplex& v) {m_CT = v;}
535 
536  private:
537  // Parameters
538 
539  // physics constant <-- not affect distributions
540  //double hbar;
541  //double GF;
542  //double Vcb;
543 
545  double m_rho12;
546 
548  double m_rhoA12;
549 
551  double m_ffR11;
552 
554  double m_ffR21;
555 
557  double m_aS1;
558 
560  double m_aR3;
561 
563  double m_mB;
564 
566  double m_mD;
567 
569  double m_mDst;
570 
572  double m_mBottom;
573 
575  double m_mCharm;
576 
578  EvtComplex m_CV1;
579 
581  EvtComplex m_CV2;
582 
584  EvtComplex m_CS1;
585 
587  EvtComplex m_CS2;
588 
590  EvtComplex m_CT;
591 
597  bool chkDhel(int Dhel) const;
598 
602  bool chkwhel(int whel) const;
603 
607  bool chktauhel(int tauhel) const;
608 
612  //bool chkcostau(double costau) const;
613 
614  };
615 
617 } // Belle 2 Namespace
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hT2
double hT2(double w) const
D* tensor form factor h_{T2}(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:442
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::dS1
double dS1(double w) const
HQET correction factor for the scalar form factor for B->Dtaunu.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:498
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::HadV1
double HadV1(int Dhel, int whel, double w) const
The function to calculate the Hadronic Amplitudes of left handed (V-A) type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:207
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::aR3
double aR3() const
HQET correction factor for the uncertainty of 1/m_Q correction.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:373
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setRhoA12
void setRhoA12(double v)
Sets the form factor parameter rho_A1^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:500
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hA1
double hA1(double w) const
HQET D* axial vector form factor h_{A1}(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:377
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::z
double z(double w) const
CLN form factor z.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:459
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hT
double hT(double w) const
D tensor form factor h_T(w) in terms of vector form factors.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:426
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_rhoA12
double m_rhoA12
Form factor slope parameters rho_A1^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:556
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::HadV2
double HadV2(int Dhel, int whel, double w) const
The function to calculate the Hadronic Amplitudes of right handed (V+A) type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:239
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helampT
double helampT(double mtau, int tauhel, int Dhel, double w, double costau) const
Helicity Amplitudes of tensor type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:347
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setCT
void setCT(const EvtComplex &v)
Sets the Wilson coeffcient CT.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:542
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::EvtBSemiTauonicHelicityAmplitudeCalculator
EvtBSemiTauonicHelicityAmplitudeCalculator()
The default constructor.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:31
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::~EvtBSemiTauonicHelicityAmplitudeCalculator
virtual ~EvtBSemiTauonicHelicityAmplitudeCalculator()
The destructor.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:74
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_mDst
double m_mDst
daughter vector (D*) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:577
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setAS1
void setAS1(double v)
Sets the form factor 1/m_Q correction parameter a_S1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:509
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_ffR11
double m_ffR11
Form factor parameter R_1(1).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:559
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_mB
double m_mB
parent (B) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:571
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setCV2
void setCV2(const EvtComplex &v)
Sets the Wilson coeffcient CV2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:533
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffS1
double ffS1(double w) const
CLN form factor S1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:471
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::wfunc
double wfunc(int Dhel, double q2) const
Calculate the velocity transfer variable w.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:432
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hT1
double hT1(double w) const
D* tensor form factor h_{T1}(w) in terms of axial vector form factors.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:433
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getRho12
double getRho12() const
parameter accessor
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:449
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_CV1
EvtComplex m_CV1
Wilson coefficient CV1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:586
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffR3
double ffR3(double w) const
CLN form factor R3.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:493
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::qh2
double qh2(int Dhel, double w) const
Function to calculate the q^2 divided by the square of parent mass (m_B^2).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:526
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::q2
double q2(int Dhel, double w) const
Function to calculate the q^2 of the decay (square of l+nu invariant mass).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:532
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::v
double v(double mtau, double q2) const
Function to calculate the tau velocity.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:520
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setMBottom
void setMBottom(double m)
Returns the bottom quark mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:524
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setMDst
void setMDst(double m)
Sets the daughter vector (D) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:521
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getMDst
double getMDst() const
Returns the daughter vector (D*) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:473
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getMD
double getMD() const
Returns the daughter scalar (D) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:470
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getMBottom
double getMBottom() const
Returns the bottom quark mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:476
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setMB
void setMB(double m)
Sets the parent (B) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:515
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hA3
double hA3(double w) const
HQET D* axial vector form factor h_{A3}(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:395
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setCV1
void setCV1(const EvtComplex &v)
Sets the Wilson coeffcient CV1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:530
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setRho12
void setRho12(double v)
Sets the form factor parameter rho_1^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:497
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getAS1
double getAS1() const
Returns form factor 1/m_Q correction factor a_S1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:461
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setR21
void setR21(double v)
Sets the form factor parameter R_2(1).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:506
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffV1
double ffV1(double w) const
CLN form factor V1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:464
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::rq
double rq() const
Ratio of the charm quark mass to the charm quark mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:390
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::r
double r(int Dhel) const
Ratio of the daughter meson mass to the parent meson.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:385
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getCS1
EvtComplex getCS1() const
Returns the Wilson coeffcient CS1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:488
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::dR3
double dR3(double w) const
HQET correction factor for the scalar form factor for B->D*taunu.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:503
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_ffR21
double m_ffR21
Form factor parameter R_2(1).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:562
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::wmax
double wmax(double mtau, int Dhel) const
Maximum value of the velocity transfer variable w.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:444
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hm
double hm(double w) const
HQET D vector form factor h_-(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:371
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getMB
double getMB() const
Returns the parent (B) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:467
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::aS1
double aS1() const
HQET correction factor for the uncertainty of 1/m_Q correction.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:368
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hV
double hV(double w) const
HQET D* axial vector form factor h_V(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:383
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getR21
double getR21() const
Returns form factor parameter R_2(1).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:458
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setAR3
void setAR3(double v)
Sets the form factor 1/m_Q correction parameter a_R3.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:512
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helampS2
double helampS2(double mtau, int tauhel, int Dhel, double w, double costau) const
Helicity Amplitudes of scalar (S-P) type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:341
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::q2max
double q2max(int Dhel) const
Maximum value of the q^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:425
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hS
double hS(double w) const
D scalar form factor h_S(w) in terms of vector form factors.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:402
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffV11
double ffV11() const
Form factor normalization factor for B->Dlnu.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:345
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_mD
double m_mD
daughter scalar (D) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:574
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::eta
double eta(int whel) const
The metric factor.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:252
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::Lep
double Lep(const double mtau, int tauhel, int whel, double q2, double costau) const
The function to calculate the Leptonic Amplitudes for B->D*taunu decay of the vector type contributio...
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:137
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getCV1
EvtComplex getCV1() const
Returns the Wilson coeffcient CV1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:482
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_aR3
double m_aR3
1/mQ correcion factor a_R3.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:568
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getR11
double getR11() const
Returns form factor parameter R_1(1).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:455
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getMCharm
double getMCharm() const
Returns the charm quark mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:479
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::wmin
double wmin() const
Minimum value of the velocity transfer variable w.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:437
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getCS2
EvtComplex getCS2() const
Returns the Wilson coeffcient CS2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:491
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::HadS2
double HadS2(int Dhel, double w) const
The function to calculate the Hadronic Amplitudes of scalar (S-P) type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:264
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setCS1
void setCS1(const EvtComplex &v)
Sets the Wilson coeffcient CS1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:536
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_aS1
double m_aS1
1/mQ correcion factor a_S1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:565
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffR1
double ffR1(double w) const
CLN form factor R1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:483
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_CS1
EvtComplex m_CS1
Wilson coefficient CS1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:592
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setCS2
void setCS2(const EvtComplex &v)
Sets the Wilson coeffcient CS2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:539
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getCV2
EvtComplex getCV2() const
Returns the Wilson coeffcient CV2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:485
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helampS1
double helampS1(double mtau, int tauhel, int Dhel, double w, double costau) const
Helicity Amplitudes of scalar (S+P) type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:335
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getCT
EvtComplex getCT() const
Returns the Wilson coeffcient CT.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:494
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::HadT
double HadT(int Dhel, int whel1, int whel2, double w) const
The function to calculate the Hadronic Amplitudes of tensor type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:274
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::q2min
double q2min(double mtau) const
Minimum value of the q^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:419
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_mBottom
double m_mBottom
b quark mass (running mass at m_b scale), used for scalar form factor term
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:580
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffA1
double ffA1(double w) const
CLN form factor A1.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:476
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::HadS1
double HadS1(int Dhel, double w) const
The function to calculate the Hadronic Amplitudes of scalar (S+P) type contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:252
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hT3
double hT3(double w) const
D* tensor form factor h_{T3}(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:450
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helampSM
double helampSM(double mtau, int tauhel, int Dhel, double w, double costau) const
Helicity Amplitudes of SM (left handed) contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:307
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setR11
void setR11(double v)
Sets the form factor parameter R_1(1).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:503
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::chkDhel
bool chkDhel(int Dhel) const
sanity checkers
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:538
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setMD
void setMD(double m)
Sets the daughter scalar (D) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:518
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helampV2
double helampV2(double mtau, int tauhel, int Dhel, double w, double costau) const
Helicity Amplitudes of right handed (V+A) contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:324
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getAR3
double getAR3() const
Returns form factor 1/m_Q correction factor a_R3.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:464
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::mD
double mD(int Dhel) const
Daughter D(*) meson mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:509
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_CS2
EvtComplex m_CS2
Wilson coefficient CS2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:595
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffA11
double ffA11() const
Form factor normalization factor for B->D*lnu.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:351
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_mCharm
double m_mCharm
c quark mass (running mass at m_b scale), used for scalar form factor term )
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:583
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helAmp
EvtComplex helAmp(double mtau, int tauhel, int Dhel, double w, double costau) const
The function calculates the helicity amplitude.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:111
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::helampV1
double helampV1(double mtau, int tauhel, int Dhel, double w, double costau) const
Helicity Amplitudes of left handed (V-A) contribution.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:318
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::getRhoA12
double getRhoA12() const
Returns form factor parameter rho_A1^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:452
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::chktauhel
bool chktauhel(int tauhel) const
Function to check if tauhel is in the valid range.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:550
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hp
double hp(double w) const
HQET D vector form factor h_+(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:364
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::chkwhel
bool chkwhel(int whel) const
Function to check if whel is in the valid range.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:544
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::setMCharm
void setMCharm(double m)
Returns the charm quark mass.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:527
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hP
double hP(double w) const
D* pseudo scalar form factor h_P(w) in terms of axial vector form factors.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:414
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::hA2
double hA2(double w) const
HQET D* axial vector form factor h_{A2}(w).
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:389
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_CT
EvtComplex m_CT
Wilson coefficient CT.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:598
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::ffR2
double ffR2(double w) const
CLN form factor R2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.cc:488
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_CV2
EvtComplex m_CV2
Wilson coefficient CV2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:589
Belle2::EvtBSemiTauonicHelicityAmplitudeCalculator::m_rho12
double m_rho12
Form factor slope parameters rho_1^2.
Definition: EvtBSemiTauonicHelicityAmplitudeCalculator.h:553