GM2Calc 2.3.0
Loading...
Searching...
No Matches
THDM.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_HPP
20#define GM2_THDM_HPP
21
23#include "SM.hpp"
24
25#include <iosfwd>
26
27#include <Eigen/Core>
28
29namespace gm2calc {
30
31namespace thdm {
32
33/**
34 * @class Config
35 * @brief Configuration options for the THDM
36 */
37struct Config {
38 bool force_output{false}; ///< force output
39 bool running_couplings{true}; ///< use running couplings
40};
41
42enum class Yukawa_type : int {
43 type_1 = 1,
44 type_2,
45 type_X,
46 type_Y,
47 aligned,
49};
50
51/** convert int to thdm::Yukawa_type */
53
56 Eigen::Matrix<double,7,1> lambda{Eigen::Matrix<double,7,1>::Zero()};
57 double tan_beta{0.0};
58 double m122{0.0};
59 double zeta_u{0.0};
60 double zeta_d{0.0};
61 double zeta_l{0.0};
62 Eigen::Matrix<double,3,3> Delta_u{Eigen::Matrix<double,3,3>::Zero()};
63 Eigen::Matrix<double,3,3> Delta_d{Eigen::Matrix<double,3,3>::Zero()};
64 Eigen::Matrix<double,3,3> Delta_l{Eigen::Matrix<double,3,3>::Zero()};
65 Eigen::Matrix<double,3,3> Pi_u{Eigen::Matrix<double,3,3>::Zero()};
66 Eigen::Matrix<double,3,3> Pi_d{Eigen::Matrix<double,3,3>::Zero()};
67 Eigen::Matrix<double,3,3> Pi_l{Eigen::Matrix<double,3,3>::Zero()};
68};
69
70struct Mass_basis {
72 double mh{0.0};
73 double mH{0.0};
74 double mA{0.0};
75 double mHp{0.0};
77 double lambda_6{0.0};
78 double lambda_7{0.0};
79 double tan_beta{0.0};
80 double m122{0.0};
81 double zeta_u{0.0};
82 double zeta_d{0.0};
83 double zeta_l{0.0};
84 Eigen::Matrix<double,3,3> Delta_u{Eigen::Matrix<double,3,3>::Zero()};
85 Eigen::Matrix<double,3,3> Delta_d{Eigen::Matrix<double,3,3>::Zero()};
86 Eigen::Matrix<double,3,3> Delta_l{Eigen::Matrix<double,3,3>::Zero()};
87 Eigen::Matrix<double,3,3> Pi_u{Eigen::Matrix<double,3,3>::Zero()};
88 Eigen::Matrix<double,3,3> Pi_d{Eigen::Matrix<double,3,3>::Zero()};
89 Eigen::Matrix<double,3,3> Pi_l{Eigen::Matrix<double,3,3>::Zero()};
90};
91
92} // namespace thdm
93
94/**
95 * @class THDM
96 * @brief Contains routines to determine the THDM parameters
97 */
98class THDM : private THDM_mass_eigenstates {
99public:
100 THDM(const thdm::Gauge_basis&, const SM& sm_ = SM{}, const thdm::Config& cfg = thdm::Config{});
101 THDM(const thdm::Mass_basis&, const SM& sm_ = SM{}, const thdm::Config& cfg = thdm::Config{});
102
103 void print(std::ostream&) const;
104
105 double get_zeta_u() const;
106 double get_zeta_d() const;
107 double get_zeta_l() const;
108
109 Eigen::Matrix<std::complex<double>,3,3> get_yuh() const;
110 Eigen::Matrix<std::complex<double>,3,3> get_yuH() const;
111 Eigen::Matrix<std::complex<double>,3,3> get_yuA() const;
112 Eigen::Matrix<std::complex<double>,3,3> get_yuHp() const;
113
114 Eigen::Matrix<std::complex<double>,3,3> get_ydh() const;
115 Eigen::Matrix<std::complex<double>,3,3> get_ydH() const;
116 Eigen::Matrix<std::complex<double>,3,3> get_ydA() const;
117 Eigen::Matrix<std::complex<double>,3,3> get_ydHp() const;
118
119 Eigen::Matrix<std::complex<double>,3,3> get_ylh() const;
120 Eigen::Matrix<std::complex<double>,3,3> get_ylH() const;
121 Eigen::Matrix<std::complex<double>,3,3> get_ylA() const;
122 Eigen::Matrix<std::complex<double>,3,3> get_ylHp() const;
123
124 const SM& get_sm() const { return sm; }
125
126 void set_tan_beta(double);
127
177
179
180private:
181 SM sm{};
183 double zeta_u{0.0}; ///< alignment parameter
184 double zeta_d{0.0}; ///< alignment parameter
185 double zeta_l{0.0}; ///< alignment parameter
186 Eigen::Matrix<double,3,3> Delta_u{Eigen::Matrix<double,3,3>::Zero()}; ///< deviation from alignment
187 Eigen::Matrix<double,3,3> Delta_d{Eigen::Matrix<double,3,3>::Zero()}; ///< deviation from alignment
188 Eigen::Matrix<double,3,3> Delta_l{Eigen::Matrix<double,3,3>::Zero()}; ///< deviation from alignment
189 thdm::Config config{}; ///< configuration options
190
191 Eigen::Matrix<double,3,1> get_mu(double) const;
192 Eigen::Matrix<double,3,1> get_md(double) const;
193 Eigen::Matrix<double,3,1> get_ml(double) const;
194 Eigen::Matrix<std::complex<double>,3,3> get_rho_u(const Eigen::Matrix<double,3,3>&) const;
195 Eigen::Matrix<std::complex<double>,3,3> get_rho_d(const Eigen::Matrix<double,3,3>&) const;
196 Eigen::Matrix<std::complex<double>,3,3> get_rho_l(const Eigen::Matrix<double,3,3>&) const;
197 void init_gauge_couplings();
198 void init_yukawas();
199 void set_basis(const thdm::Gauge_basis&);
200 void set_basis(const thdm::Mass_basis&);
201 void validate() const;
202 const char* yukawa_type_to_string() const;
203};
204
205/// streaming operator
206std::ostream& operator<<(std::ostream&, const THDM&);
207
208} // namespace gm2calc
209
210#endif
contains class for general THDM model with routines needed to solve EWSB and determine the masses and...
model class with routines for determing masses and mixinga and EWSB
const THDM_problems & get_problems() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Ue() const
const Eigen::Array< double, 3, 1 > & get_MFe() const
charged lepton masses
const Eigen::Array< double, 2, 1 > & get_MHm() const
Goldstone and charged Higgs boson masses (in that order)
double get_LambdaFive() const
capital Lambda5, Eq (14) arxiv:1607.06292
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Vd() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Ud() const
double get_sin_beta_minus_alpha() const
sin(beta - alpha_h)
const Eigen::Matrix< double, 2, 2 > & get_ZH() const
CP-even Higgs boson mixing matrix.
double get_MVZ() const
Z boson mass.
double get_cos_beta_minus_alpha() const
cos(beta - alpha_h)
const Eigen::Array< double, 3, 1 > & get_MFu() const
up-type quark masses
double get_tan_beta() const
tan(beta) = ratio of VEVs v2/v1
const Eigen::Array< double, 3, 1 > & get_MFv() const
neutrino masses
const Eigen::Matrix< double, 2, 2 > & get_ZA() const
CP-odd Higgs boson mixing matrix.
double get_alpha_h() const
CP-even Higgs mixing angle.
double get_v_sqr() const
squared SM-like VEV
double get_beta() const
CP-odd and charged Higgs mixing angle.
const Eigen::Array< double, 2, 1 > & get_MAh() const
Goldstone and CP-odd Higgs boson masses (in that order)
const Eigen::Array< double, 3, 1 > & get_MFd() const
down-type quark masses
double get_alpha_em() const
electromagnetic coupling
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Ve() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Vu() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Uu() const
const Eigen::Array< double, 2, 1 > & get_Mhh() const
CP-even Higgs boson masses.
double get_LambdaSixSeven() const
(Lambda_{567} - Lambda_{5})(tan(b) - 1/tan(b))
const Eigen::Matrix< double, 2, 2 > & get_ZP() const
charged Higgs boson mixing matrix
double get_eta() const
deviation of CP-even Higgs mixing angle from SM limit
double get_MVWm() const
W boson mass.
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Pi_d() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Gamma_d() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Pi_l() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Pi_u() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Gamma_l() const
const Eigen::Matrix< std::complex< double >, 3, 3 > & get_Gamma_u() const
Contains routines to determine the THDM parameters.
Definition THDM.hpp:98
Eigen::Matrix< std::complex< double >, 3, 3 > get_ydH() const
Definition THDM.cpp:334
Eigen::Matrix< std::complex< double >, 3, 3 > get_ylH() const
Definition THDM.cpp:366
double get_zeta_u() const
Table 1, arxiv:1607.06292.
Definition THDM.cpp:200
double get_zeta_l() const
Table 1, arxiv:1607.06292.
Definition THDM.cpp:237
Eigen::Matrix< std::complex< double >, 3, 3 > get_yuA() const
Definition THDM.cpp:312
Eigen::Matrix< std::complex< double >, 3, 3 > get_ydA() const
Definition THDM.cpp:344
Eigen::Matrix< std::complex< double >, 3, 3 > get_ylHp() const
Definition THDM.cpp:382
Eigen::Matrix< std::complex< double >, 3, 3 > get_yuH() const
Definition THDM.cpp:302
const SM & get_sm() const
Definition THDM.hpp:124
void print(std::ostream &) const
Definition THDM.cpp:546
Eigen::Matrix< std::complex< double >, 3, 3 > get_yuHp() const
Definition THDM.cpp:318
Eigen::Matrix< std::complex< double >, 3, 3 > get_yuh() const
Definition THDM.cpp:292
Eigen::Matrix< std::complex< double >, 3, 3 > get_ylh() const
Definition THDM.cpp:356
Eigen::Matrix< std::complex< double >, 3, 3 > get_ydHp() const
Definition THDM.cpp:350
Eigen::Matrix< std::complex< double >, 3, 3 > get_ydh() const
Definition THDM.cpp:324
Eigen::Matrix< std::complex< double >, 3, 3 > get_ylA() const
Definition THDM.cpp:376
double get_zeta_d() const
Table 1, arxiv:1607.06292.
Definition THDM.cpp:217
void set_tan_beta(double)
Definition THDM.cpp:541
Yukawa_type int_to_cpp_yukawa_type(int yukawa_type)
convert int to thdm::Yukawa_type
Definition THDM.cpp:45
std::ostream & operator<<(std::ostream &os, const MSSMNoFV_onshell &model)
streaming operator
Configuration options for the THDM.
Definition THDM.hpp:37
bool running_couplings
use running couplings
Definition THDM.hpp:39
bool force_output
force output
Definition THDM.hpp:38
Eigen::Matrix< double, 3, 3 > Pi_u
Definition THDM.hpp:65
Eigen::Matrix< double, 3, 3 > Delta_u
Definition THDM.hpp:62
Eigen::Matrix< double, 3, 3 > Pi_l
Definition THDM.hpp:67
Eigen::Matrix< double, 3, 3 > Delta_l
Definition THDM.hpp:64
Eigen::Matrix< double, 7, 1 > lambda
Definition THDM.hpp:56
Eigen::Matrix< double, 3, 3 > Pi_d
Definition THDM.hpp:66
Yukawa_type yukawa_type
Definition THDM.hpp:55
Eigen::Matrix< double, 3, 3 > Delta_d
Definition THDM.hpp:63
Yukawa_type yukawa_type
Definition THDM.hpp:71
Eigen::Matrix< double, 3, 3 > Delta_u
Definition THDM.hpp:84
Eigen::Matrix< double, 3, 3 > Pi_u
Definition THDM.hpp:87
Eigen::Matrix< double, 3, 3 > Delta_d
Definition THDM.hpp:85
Eigen::Matrix< double, 3, 3 > Pi_d
Definition THDM.hpp:88
Eigen::Matrix< double, 3, 3 > Pi_l
Definition THDM.hpp:89
Eigen::Matrix< double, 3, 3 > Delta_l
Definition THDM.hpp:86