flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
library_softsusy.cpp
Go to the documentation of this file.
1// ====================================================================
2// This file is part of FlexibleSUSY.
3//
4// FlexibleSUSY 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// FlexibleSUSY 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 FlexibleSUSY. If not, see
16// <http://www.gnu.org/licenses/>.
17// ====================================================================
18
19#include "library_softsusy.hpp"
20#include "numerics.h"
21#include <limits>
22
23#define NAN_Q std::numeric_limits<double>::quiet_NaN()
24#define UNDEFINED(R, ARGS, NAME) \
25 std::complex<double> Softsusy::NAME ARGS noexcept { return {NAN_Q, NAN_Q}; }
26
27namespace flexiblesusy
28{
29namespace looplibrary
30{
31
32std::complex<double> Softsusy::A0(A_ARGS) noexcept
33{
34 double m2 = m02_in.real();
35 double q2 = scl2_in;
36
37 return {softsusy::rea0(m2, q2), 0.0};
38}
39
40std::complex<double> Softsusy::B0(B_ARGS) noexcept
41{
42 double p = std::sqrt(p10_in.real());
43 double m1 = std::sqrt(m02_in.real());
44 double m2 = std::sqrt(m12_in.real());
45 double q = std::sqrt(scl2_in);
46
47 return {softsusy::b0(p, m1, m2, q), 0.0};
48}
49
50std::complex<double> Softsusy::B1(B_ARGS) noexcept
51{
52 double p = std::sqrt(p10_in.real());
53 double m1 = std::sqrt(m02_in.real());
54 double m2 = std::sqrt(m12_in.real());
55 double q = std::sqrt(scl2_in);
56
57 return {(-1) * softsusy::b1(p, m1, m2, q), 0.0};
58}
59
60std::complex<double> Softsusy::B00(B_ARGS) noexcept
61{
62 double p = std::sqrt(p10_in.real());
63 double m1 = std::sqrt(m02_in.real());
64 double m2 = std::sqrt(m12_in.real());
65 double q = std::sqrt(scl2_in);
66
67 return {softsusy::b22(p, m1, m2, q), 0.0};
68}
69
70std::complex<double> Softsusy::C0(C_ARGS) noexcept
71{
72 double m1 = std::sqrt(m02_in.real());
73 double m2 = std::sqrt(m12_in.real());
74 double m3 = std::sqrt(m22_in.real());
75
76 return {softsusy::c0(m1, m2, m3), 0.0};
77}
78
79std::complex<double> Softsusy::C00(C_ARGS) noexcept
80{
81 double m1 = std::sqrt(m02_in.real());
82 double m2 = std::sqrt(m12_in.real());
83 double m3 = std::sqrt(m22_in.real());
84 double q = std::sqrt(scl2_in);
85
86 return {softsusy::c00(m1, m2, m3, q), 0.0};
87}
88
89BOOST_PP_SEQ_FOR_EACH(UNDEFINED, (C_ARGS), (C1)(C2)(C11)(C12)(C22))
90
91std::complex<double> Softsusy::D0(D_ARGS) noexcept
92{
93 double m1 = std::sqrt(m02_in.real());
94 double m2 = std::sqrt(m12_in.real());
95 double m3 = std::sqrt(m22_in.real());
96 double m4 = std::sqrt(m32_in.real());
97
98 return {softsusy::d0(m1, m2, m3, m4), 0.0};
99}
100
101std::complex<double> Softsusy::D00(D_ARGS) noexcept
102{
103 double m1 = std::sqrt(m02_in.real());
104 double m2 = std::sqrt(m12_in.real());
105 double m3 = std::sqrt(m22_in.real());
106 double m4 = std::sqrt(m32_in.real());
107
108 return {softsusy::d27(m1, m2, m3, m4), 0.0};
109}
110
111BOOST_PP_SEQ_FOR_EACH(UNDEFINED, (D_ARGS),
112 (D1)(D11)(D12)(D13)(D2)(D22)(D23)(D3)(D33))
113
114void Softsusy::A(Acoeff_t& a, A_ARGS) noexcept
115{
116 double m = std::sqrt(m02_in.real());
117 double q = std::sqrt(scl2_in);
118
119 a.at(0) = {softsusy::a0(m, q), 0.0};
120}
121
122void Softsusy::B(Bcoeff_t& b, B_ARGS) noexcept
123{
124 double p = std::sqrt(p10_in.real());
125 double m1 = std::sqrt(m02_in.real());
126 double m2 = std::sqrt(m12_in.real());
127 double q = std::sqrt(scl2_in);
128
129 b.at(0) = {softsusy::b0(p, m1, m2, q), 0.0};
130 b.at(1) = {(-1) * softsusy::b1(p, m1, m2, q), 0.0};
131 b.at(2) = {softsusy::b22(p, m1, m2, q), 0.0};
132}
133
134void Softsusy::C(Ccoeff_t& c, C_ARGS) noexcept
135{
136 double m1 = std::sqrt(m02_in.real());
137 double m2 = std::sqrt(m12_in.real());
138 double m3 = std::sqrt(m22_in.real());
139 double q = std::sqrt(scl2_in);
140 std::complex<double> undefined = {NAN_Q, NAN_Q};
141
142 c.at(0) = {softsusy::c0(m1, m2, m3), 0.0};
143 c.at(1) = undefined;
144 c.at(2) = undefined;
145 c.at(3) = {softsusy::c00(m1, m2, m3, q), 0.0};
146 c.at(4) = undefined;
147 c.at(5) = undefined;
148 c.at(6) = undefined;
149}
150
151void Softsusy::D(Dcoeff_t& d, D_ARGS) noexcept
152{
153 double m1 = std::sqrt(m02_in.real());
154 double m2 = std::sqrt(m12_in.real());
155 double m3 = std::sqrt(m22_in.real());
156 double m4 = std::sqrt(m32_in.real());
157 std::complex<double> undefined = {NAN_Q, NAN_Q};
158
159 d.at(0) = {softsusy::d0(m1, m2, m3, m4), 0.0};
160 d.at(1) = undefined;
161 d.at(2) = undefined;
162 d.at(3) = undefined;
163 d.at(4) = {softsusy::d27(m1, m2, m3, m4), 0.0};
164 d.at(5) = undefined;
165 d.at(6) = undefined;
166 d.at(7) = undefined;
167 d.at(8) = undefined;
168 d.at(9) = undefined;
169 d.at(10) = undefined;
170}
171
172} // namespace looplibrary
173} // namespace flexiblesusy
void B(Bcoeff_t &, BOOST_PP_SEQ_FOR_EACH(ARGS_TYPE,,(p10_in)(m02_in)(m12_in)) double scl2_in) noexcept override
void C(Ccoeff_t &, BOOST_PP_SEQ_FOR_EACH(ARGS_TYPE,,(p10_in)(p21_in)(p20_in)(m02_in)(m12_in)(m22_in)) double scl2_in) noexcept override
void D(Dcoeff_t &, BOOST_PP_SEQ_FOR_EACH(ARGS_TYPE,,(p10_in)(p21_in)(p32_in)(p30_in)(p20_in)(p31_in)(m02_in)(m12_in)(m22_in)(m32_in)) double scl2_in) noexcept override
void A(Acoeff_t &, BOOST_PP_SEQ_FOR_EACH(ARGS_TYPE,,(m02_in)) double scl2_in) noexcept override
#define NAN_Q
#define UNDEFINED(R, ARGS, NAME)
#define B_ARGS
#define C_ARGS
#define D_ARGS
#define A_ARGS
std::array< std::complex< double >, BOOST_PP_SEQ_SIZE((0)(1)(2)(00)(11)(12)(22)) > Ccoeff_t
std::array< std::complex< double >, BOOST_PP_SEQ_SIZE((0)) > Acoeff_t
std::array< std::complex< double >, BOOST_PP_SEQ_SIZE((0)(1)(2)(3)(00)(11)(12)(13)(22)(23)(33)) > Dcoeff_t
std::array< std::complex< double >, BOOST_PP_SEQ_SIZE((0)(1)(00)) > Bcoeff_t
double C0(double m1, double m2, double m3) noexcept
(arguments are interpreted as unsquared)
double D0(double m1, double m2, double m3, double m4) noexcept
(arguments are interpreted as unsquared)
double B0(double m1, double m2, double scale) noexcept
(arguments are interpreted as unsquared)
double b1(double p, double m1, double m2, double q) noexcept
Note that b1 is NOT symmetric in m1 <-> m2!!!
Definition: numerics.cpp:223
double c00(double m1, double m2, double m3, double q) noexcept
Definition: numerics.cpp:458
double rea0(double x, double q) noexcept
Definition: numerics.cpp:139
double c0(double m1, double m2, double m3) noexcept
Definition: numerics.cpp:373
double a0(double m, double q) noexcept
Definition: numerics.cpp:125
double b22(double p, double m1, double m2, double q) noexcept
Definition: numerics.cpp:275
double d27(double m1, double m2, double m3, double m4) noexcept
Definition: numerics.cpp:362
double b0(double p, double m1, double m2, double q) noexcept
Definition: numerics.cpp:174
double d0(double m1, double m2, double m3, double m4) noexcept
Definition: numerics.cpp:322
loop functions