GM2Calc 2.3.0
Loading...
Searching...
No Matches
gm2_2loop_helpers.hpp
Go to the documentation of this file.
1// ====================================================================
2// This file is part of GM2Calc.
3//
4// GM2Calc is free software: you can redistribute it and/or modify
5// it under the terms of the GNU General Public License as published
6// by the Free Software Foundation, either version 3 of the License,
7// or (at your option) any later version.
8//
9// GM2Calc is distributed in the hope that it will be useful, but
10// WITHOUT ANY WARRANTY; without even the implied warranty of
11// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12// General Public License for more details.
13//
14// You should have received a copy of the GNU General Public License
15// along with GM2Calc. If not, see
16// <http://www.gnu.org/licenses/>.
17// ====================================================================
18
19#ifndef GM2_THDM_2LOOP_HELPERS_HPP
20#define GM2_THDM_2LOOP_HELPERS_HPP
21
22#include <Eigen/Core>
23
24namespace gm2calc {
25
26namespace thdm {
27
28/// parameters to be passed to the bosonic contribution functions
30 double alpha_em{};///< electromagnetic coupling
31 double mm{}; ///< muon mass for prefactor
32 double mw{}; ///< W boson mass
33 double mz{}; ///< Z boson mass
34 double mhSM{}; ///< SM Higgs boson mass
35 double mA{}; ///< CP-odd Higgs boson mass
36 double mHp{}; ///< charged Higgs boson mass
37 Eigen::Matrix<double,2,1> mh{Eigen::Matrix<double,2,1>::Zero()}; ///< CP-even Higgs bosons mass
38 double tb{}; ///< tan(beta)
39 double zetal{}; ///< zeta_l
40 double cos_beta_minus_alpha{}; ///< cos(beta - alpha_h)
41 double lambda5{}; ///< Lambda_5
42 double lambda67{}; ///< difference (Lambda_567 - Lambda_5)
43};
44
45/// parameters to be passed to the fermionic contribution functions
47 double alpha_em{}; ///< electromagnetic coupling
48 double mm{}; ///< muon mass for prefactor
49 double mw{}; ///< W boson mass
50 double mz{}; ///< Z boson mass
51 double mhSM{}; ///< SM Higgs boson mass
52 double mA{}; ///< CP-odd Higgs boson mass
53 double mHp{}; ///< charged Higgs boson mass
54 Eigen::Matrix<double,2,1> mh{Eigen::Matrix<double,2,1>::Zero()}; ///< CP-even Higgs bosons mass
55 Eigen::Matrix<double,3,1> ml{Eigen::Matrix<double,3,1>::Zero()}; ///< down-type lepton masses
56 Eigen::Matrix<double,3,1> mu{Eigen::Matrix<double,3,1>::Zero()}; ///< up-type quark masses
57 Eigen::Matrix<double,3,1> md{Eigen::Matrix<double,3,1>::Zero()}; ///< down-type quark masses
58 Eigen::Matrix<std::complex<double>,3,3> yuh{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={u,c,t} and S=h
59 Eigen::Matrix<std::complex<double>,3,3> yuH{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={u,c,t} and S=H
60 Eigen::Matrix<std::complex<double>,3,3> yuA{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={u,c,t} and S=A
61 Eigen::Matrix<std::complex<double>,3,3> yuHp{Eigen::Matrix<std::complex<double>,3,3>::Zero()};///< y_f^S coefficients with f={u,c,t} and S=H^+
62 Eigen::Matrix<std::complex<double>,3,3> ydh{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={d,s,b} and S=h
63 Eigen::Matrix<std::complex<double>,3,3> ydH{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={d,s,b} and S=H
64 Eigen::Matrix<std::complex<double>,3,3> ydA{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={d,s,b} and S=A
65 Eigen::Matrix<std::complex<double>,3,3> ydHp{Eigen::Matrix<std::complex<double>,3,3>::Zero()};///< y_f^S coefficients with f={d,s,b} and S=H^+
66 Eigen::Matrix<std::complex<double>,3,3> ylh{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={e,m,τ} and S=h
67 Eigen::Matrix<std::complex<double>,3,3> ylH{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={e,m,τ} and S=H
68 Eigen::Matrix<std::complex<double>,3,3> ylA{Eigen::Matrix<std::complex<double>,3,3>::Zero()}; ///< y_f^S coefficients with f={e,m,τ} and S=A
69 Eigen::Matrix<std::complex<double>,3,3> ylHp{Eigen::Matrix<std::complex<double>,3,3>::Zero()};///< y_f^S coefficients with f={e,m,τ} and S=H^+
70 Eigen::Matrix<std::complex<double>,3,3> vckm{Eigen::Matrix<std::complex<double>,3,3>::Identity()};///< CKM matrix
71};
72
73// === 2-loop bosonic contributions ===
74
75double amu2L_B(const THDM_B_parameters&) noexcept;
76
77// routines for sub-expressions
78
79double amu2L_B_EWadd(const THDM_B_parameters&) noexcept;
80double amu2L_B_nonYuk(const THDM_B_parameters&) noexcept;
81double amu2L_B_Yuk(const THDM_B_parameters&) noexcept;
82
83// === 2-loop fermionic contributions ===
84
85double amu2L_F(const THDM_F_parameters&) noexcept;
86
87// routines for sub-expressions
88
89double amu2L_F_charged(const THDM_F_parameters&) noexcept;
90double amu2L_F_neutral(const THDM_F_parameters&) noexcept;
91
92/// Eq (53), arxiv:1607.06292, f = u, S = h or H
93double fuS(double ms2, double mu2, double mw2, double mz2) noexcept;
94/// Eq (53), arxiv:1607.06292, f = d, S = h or H
95double fdS(double ms2, double md2, double mw2, double mz2) noexcept;
96/// Eq (53), arxiv:1607.06292, f = l, S = h or H
97double flS(double ms2, double ml2, double mw2, double mz2) noexcept;
98/// Eq (53), arxiv:1607.06292, f = u, S = A
99double fuA(double ms2, double mu2, double mw2, double mz2) noexcept;
100/// Eq (53), arxiv:1607.06292, f = d, S = A
101double fdA(double ms2, double md2, double mw2, double mz2) noexcept;
102/// Eq (53), arxiv:1607.06292, f = l, S = A
103double flA(double ms2, double ml2, double mw2, double mz2) noexcept;
104/// Eq (59), arxiv:1607.06292, S = H^\pm, f = u
105double fuHp(double ms2, double md2, double mu2, double mw2, double mz2) noexcept;
106/// Eq (59), arxiv:1607.06292, S = H^\pm, f = d
107double fdHp(double ms2, double md2, double mu2, double mw2, double mz2) noexcept;
108/// Eq (59), arxiv:1607.06292, S = H^\pm, f = l
109double flHp(double ms2, double ml2, double mw2, double mz2) noexcept;
110
111} // namespace thdm
112
113} // namespace gm2calc
114
115#endif
double mw2
squared W boson mass
double mz2
squared Z boson mass
double fdA(double ms2, double md2, double mw2, double mz2) noexcept
Eq (53), arxiv:1607.06292, f = d, S = A.
double flHp(double ms2, double ml2, double mw2, double mz2) noexcept
Eq (59), arxiv:1607.06292, S = H^\pm, f = l.
double flA(double ms2, double ml2, double mw2, double mz2) noexcept
Eq (53), arxiv:1607.06292, f = l, S = A.
double fuHp(double ms2, double md2, double mu2, double mw2, double mz2) noexcept
Eq (59), arxiv:1607.06292, S = H^\pm, f = u.
double amu2L_B(const THDM_B_parameters &thdm) noexcept
Calculates the sum of the 2-loop bosonic contributions.
double amu2L_F(const THDM_F_parameters &thdm) noexcept
Calculates the sum of the 2-loop fermionic contributions with neutral and charged Higgs bosons.
double fdHp(double ms2, double md2, double mu2, double mw2, double mz2) noexcept
Eq (59), arxiv:1607.06292, S = H^\pm, f = d.
double fdS(double ms2, double md2, double mw2, double mz2) noexcept
Eq (53), arxiv:1607.06292, f = d, S = h or H.
double fuA(double ms2, double mu2, double mw2, double mz2) noexcept
Eq (53), arxiv:1607.06292, f = u, S = A.
double flS(double ms2, double ml2, double mw2, double mz2) noexcept
Eq (53), arxiv:1607.06292, f = l, S = h or H.
double amu2L_B_EWadd(const THDM_B_parameters &thdm) noexcept
Calculates 2-loop bosonic pure electroweak contributions.
double fuS(double ms2, double mu2, double mw2, double mz2) noexcept
Eq (53), arxiv:1607.06292, f = u, S = h or H.
double amu2L_F_neutral(const THDM_F_parameters &thdm) noexcept
Calculates 2-loop fermionic contributions with neutral Higgs bosons.
double amu2L_B_nonYuk(const THDM_B_parameters &thdm) noexcept
Calculates 2-loop bosonic non-Yukawa contributions.
double amu2L_B_Yuk(const THDM_B_parameters &thdm) noexcept
Calculates 2-loop bosonic Yukawa contributions.
double amu2L_F_charged(const THDM_F_parameters &thdm) noexcept
Calculates 2-loop fermionic contributions with charged Higgs boson.
void svd_eigen(const Eigen::Matrix< Scalar, M, N > &m, Eigen::Array< Real,(((M)<(N)) ?(M) :(N)), 1 > &s, Eigen::Matrix< Scalar, M, M > *u, Eigen::Matrix< Scalar, N, N > *vh)
parameters to be passed to the bosonic contribution functions
double mHp
charged Higgs boson mass
Eigen::Matrix< double, 2, 1 > mh
CP-even Higgs bosons mass.
double alpha_em
electromagnetic coupling
double mm
muon mass for prefactor
double mA
CP-odd Higgs boson mass.
double mhSM
SM Higgs boson mass.
double lambda67
difference (Lambda_567 - Lambda_5)
double cos_beta_minus_alpha
cos(beta - alpha_h)
parameters to be passed to the fermionic contribution functions
Eigen::Matrix< std::complex< double >, 3, 3 > ylH
y_f^S coefficients with f={e,m,τ} and S=H
Eigen::Matrix< std::complex< double >, 3, 3 > ylA
y_f^S coefficients with f={e,m,τ} and S=A
double alpha_em
electromagnetic coupling
Eigen::Matrix< std::complex< double >, 3, 3 > ydHp
y_f^S coefficients with f={d,s,b} and S=H^+
Eigen::Matrix< std::complex< double >, 3, 3 > vckm
CKM matrix.
Eigen::Matrix< std::complex< double >, 3, 3 > ydh
y_f^S coefficients with f={d,s,b} and S=h
Eigen::Matrix< std::complex< double >, 3, 3 > ylh
y_f^S coefficients with f={e,m,τ} and S=h
double mA
CP-odd Higgs boson mass.
Eigen::Matrix< std::complex< double >, 3, 3 > yuH
y_f^S coefficients with f={u,c,t} and S=H
Eigen::Matrix< double, 2, 1 > mh
CP-even Higgs bosons mass.
Eigen::Matrix< double, 3, 1 > ml
down-type lepton masses
Eigen::Matrix< std::complex< double >, 3, 3 > yuh
y_f^S coefficients with f={u,c,t} and S=h
Eigen::Matrix< double, 3, 1 > mu
up-type quark masses
Eigen::Matrix< double, 3, 1 > md
down-type quark masses
Eigen::Matrix< std::complex< double >, 3, 3 > ylHp
y_f^S coefficients with f={e,m,τ} and S=H^+
Eigen::Matrix< std::complex< double >, 3, 3 > ydA
y_f^S coefficients with f={d,s,b} and S=A
double mHp
charged Higgs boson mass
Eigen::Matrix< std::complex< double >, 3, 3 > ydH
y_f^S coefficients with f={d,s,b} and S=H
double mhSM
SM Higgs boson mass.
double mm
muon mass for prefactor
Eigen::Matrix< std::complex< double >, 3, 3 > yuHp
y_f^S coefficients with f={u,c,t} and S=H^+
Eigen::Matrix< std::complex< double >, 3, 3 > yuA
y_f^S coefficients with f={u,c,t} and S=A