GM2Calc 2.3.0
Loading...
Searching...
No Matches
example-gm2calc.cpp
Go to the documentation of this file.
6
7#include <cstdio>
8#include <cmath>
9
12
13 const double Pi = 3.141592653589793;
14 const Eigen::Matrix<double,3,3> UnitMatrix
15 = Eigen::Matrix<double,3,3>::Identity();
16
17 // fill SM parameters
18 model.set_alpha_MZ(0.0077552); // 1L
19 model.set_alpha_thompson(0.00729735); // 2L
20 model.set_g3(std::sqrt(4 * Pi * 0.1184)); // 2L
21 model.get_physical().MFt = 173.34; // 2L
22 model.get_physical().MFb = 4.18; // 2L, mb(mb) MS-bar
23 model.get_physical().MFm = 0.1056583715; // 1L
24 model.get_physical().MFtau = 1.777; // 2L
25 model.get_physical().MVWm = 80.385; // 1L
26 model.get_physical().MVZ = 91.1876; // 1L
27
28 // fill DR-bar parameters
29 model.set_TB(10); // 1L
30 model.set_Ae(1,1,0); // 1L
31
32 // fill on-shell parameters
33 model.set_Mu(350); // 1L
34 model.set_MassB(150); // 1L
35 model.set_MassWB(300); // 1L
36 model.set_MassG(1000); // 2L
37 model.set_mq2(500 * 500 * UnitMatrix); // 2L
38 model.set_ml2(500 * 500 * UnitMatrix); // 1L(smuon)/2L
39 model.set_md2(500 * 500 * UnitMatrix); // 2L
40 model.set_mu2(500 * 500 * UnitMatrix); // 2L
41 model.set_me2(500 * 500 * UnitMatrix); // 1L(smuon)/2L
42 model.set_Au(2,2,0); // 2L
43 model.set_Ad(2,2,0); // 2L
44 model.set_Ae(2,2,0); // 2L
45 model.set_MA0(1500); // 2L
46 model.set_scale(454.7); // 2L
47
48 // calculate mass spectrum
49 model.calculate_masses();
50
51 return model;
52}
53
54int main() {
55 try {
57
58 const double amu =
61
62 const double delta_amu =
64
65 std::printf("amu = %.5e +- %.5e\n", amu, delta_amu);
66 } catch (const gm2calc::Error& e) {
67 std::printf("%s\n", e.what());
68 }
69
70 return 0;
71}
const MSSMNoFV_onshell_physical & get_physical() const
void set_md2(const Eigen::Matrix< double, 3, 3 > &md2_)
void set_mu2(const Eigen::Matrix< double, 3, 3 > &mu2_)
void set_mq2(const Eigen::Matrix< double, 3, 3 > &mq2_)
void set_ml2(const Eigen::Matrix< double, 3, 3 > &ml2_)
void set_me2(const Eigen::Matrix< double, 3, 3 > &me2_)
contains the MSSMNoFV parameters in the on-shell scheme
void calculate_masses()
calculate SUSY masses in mixed on-shell/DR-bar scheme from given on-shell/DR-bar parameters
void set_Au(const Eigen::Matrix< double, 3, 3 > &A)
soft-breaking trilinear up-type squark coupling
void set_TB(double)
set tan(beta)
void set_alpha_MZ(double)
set alpha(MZ) w/o hadronic corrections
void set_MA0(double m)
set CP-odd Higgs pole mass
void set_Ad(const Eigen::Matrix< double, 3, 3 > &A)
soft-breaking trilinear down-type squark coupling
void set_alpha_thompson(double)
set alpha in the Thomson limit
void set_Ae(const Eigen::Matrix< double, 3, 3 > &A)
soft-breaking trilinear on-shell down-type slepton coupling
gm2calc::MSSMNoFV_onshell setup()
int main()
double calculate_amu_2loop(const MSSMNoFV_onshell &model)
Calculates best 2-loop SUSY contribution to a_mu with tan(beta) resummation.
double calculate_uncertainty_amu_2loop(const THDM &, double, double)
calculates uncertainty for amu(2-loop)
double calculate_amu_1loop(const MSSMNoFV_onshell &model)
Calculates full 1-loop SUSY contribution to (g-2), Eq (45) of arXiv:hep-ph/0609168.
Definition gm2_1loop.cpp:74