a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
GeneralTHDMZ2Runner.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2016 HEPfit Collaboration
3 *
4 *
5 * For the licensing terms see doc/COPYING.
6 */
7
9#include "GeneralTHDMcache.h"
10#include <gsl/gsl_errno.h>
11#include <gsl/gsl_matrix.h>
12#include <gsl/gsl_odeiv2.h>
13
15: myGTHDMZ2(static_cast<const GeneralTHDMZ2*> (&SM_i)), myZ2_at_Q(3, 5, 0.)
16{}
17
19{};
20
21int RGEsZ2(double t, const double y[], double beta[], void *flags)
22{
23 (void)(t); /* avoid unused parameter warning */
24 int flag = *(int *)flags;
25 int ord = flag%3;
26 int type = flag-ord;
27
28 double Yb1 = 0;
29 double Yb2 = 0;
30 double Ytau1 = 0;
31 double Ytau2 = 0;
32 double g1 = y[0];
33 double g2 = y[1];
34 double g3 = y[2];
35 double Yt = y[3];
36
37 if( type == 0 ) {//type I
38 Yb2 = y[4];
39 Ytau2 = y[5];
40 }
41 else if( type == 3 ) {//type II
42 Yb1 = y[4];
43 Ytau1 = y[5];
44 }
45 else if( type == 6 ) {//type X
46 Yb2 = y[4];
47 Ytau1 = y[5];
48 }
49 else if( type == 9 ) {//type Y
50 Yb1 = y[4];
51 Ytau2 = y[5];
52 }
53 else if( type == 12 ) {//inert
54 Yt = 0.;
55 }
56 else {
57 throw std::runtime_error("RGEsZ2: type should only be any of 0, 3, 6, 9, 12");
58 }
59
60 double m11_2=y[6];
61 double m22_2=y[7];
62 double m12_2=y[8];
63 double la1=y[9];
64 double la2=y[10];
65 double la3=y[11];
66 double la4=y[12];
67 double la5=y[13];
68
69 double pi = M_PI;
70
71 //beta_g1
72 beta[0] = 7.0*g1*g1*g1/(16.0*pi*pi);
73 //beta_g2
74 beta[1] = -3.0*g2*g2*g2/(16.0*pi*pi);
75 //beta_g3
76 beta[2] = -7.0*g3*g3*g3/(16.0*pi*pi);
77 //beta_Yt
78 beta[3] = Yt*(-17.0*g1*g1+3.0*(-9.0*g2*g2-32.0*g3*g3+2.0*Yb1*Yb1+6.0*Yb2*Yb2+18.0*Yt*Yt+4.0*Ytau2*Ytau2))/(192.0*pi*pi);
79 //beta_Yb
80 beta[4] = ((-5.0*g1*g1-27.0*g2*g2-96.0*g3*g3)*(Yb1+Yb2)+(6.0*Yb1+18.0*Yb2)*Yt*Yt
81 +54.0*(Yb1*Yb1*Yb1+Yb2*Yb2*Yb2+Yb1*Yb2*Yb2)+12.0*(Yb1*Ytau1*Ytau1+Yb2*Ytau2*Ytau2))/(192.0*pi*pi);
82 //beta_Ytau
83 beta[5] = ((-15.0*g1*g1-9.0*g2*g2)*(Ytau1+Ytau2)+Ytau2*(12.0*Yb2*Yb2+12.0*Yt*Yt+10.0*Ytau2*Ytau2)
84 +Ytau1*(12.0*Yb1*Yb1+10.0*Ytau1*Ytau1))/(64.0*pi*pi);
85 //beta_m11_2
86 beta[6] = (-3.0*g1*g1*m11_2-9.0*g2*g2*m11_2
87 +4.0*(m11_2*(3.0*Yb1*Yb1+Ytau1*Ytau1+3.0*la1)+m22_2*(2.0*la3+la4)))/(32.0*pi*pi);
88 //beta_m22_2
89 beta[7] = (-3.0*g1*g1*m22_2-9.0*g2*g2*m22_2
90 +4.0*(m22_2*(3.0*Yb2*Yb2+3.0*Yt*Yt+Ytau2*Ytau2+3.0*la2)+m11_2*(2.0*la3+la4)))/(32.0*pi*pi);
91 //beta_m12_2
92 beta[8] = m12_2*(-3.0*g1*g1-9.0*g2*g2
93 +2.0*(3.0*Yb1*Yb1+3.0*Yb2*Yb2+3.0*Yt*Yt+Ytau1*Ytau1+Ytau2*Ytau2+2.0*la3+4.0*la4+6.0*la5))/(32.0*pi*pi);
94 //beta_lambda_1
95 beta[9] = (3.0*g1*g1*g1*g1+9.0*g2*g2*g2*g2+6.0*g1*g1*(g2*g2-2.0*la1)-36.0*g2*g2*la1
96 +8.0*(6.0*la1*la1+2.0*la3*la3+2.0*la3*la4+la4*la4+la5*la5)
97 -16.0*(3.0*Yb1*Yb1*Yb1*Yb1+Ytau1*Ytau1*Ytau1*Ytau1)+16.0*la1*(3.0*Yb1*Yb1+Ytau1*Ytau1))/(64.0*pi*pi);
98 //beta_lambda_2
99 beta[10] = (3.0*g1*g1*g1*g1+9.0*g2*g2*g2*g2+6.0*g1*g1*(g2*g2-2.0*la2)-36.0*g2*g2*la2
100 -8.0*(6.0*Yb2*Yb2*Yb2*Yb2+6.0*Yt*Yt*Yt*Yt+2.0*Ytau2*Ytau2*Ytau2*Ytau2-6.0*Yb2*Yb2*la2
101 -6.0*Yt*Yt*la2-2.0*Ytau2*Ytau2*la2-6.0*la2*la2-2.0*la3*la3-2.0*la3*la4-la4*la4-la5*la5))/(64.0*pi*pi);
102 //beta_lambda_3
103 beta[11] = (3.0*g1*g1*g1*g1+9.0*g2*g2*g2*g2-36.0*g2*g2*la3-6.0*g1*g1*(g2*g2+2.0*la3)
104 +8.0*(3.0*Yb2*Yb2*la3+3.0*Yt*Yt*la3+Ytau1*Ytau1*la3+Ytau2*Ytau2*la3+3.0*la1*la3+3.0*la2*la3+2.0*la3*la3
105 +Yb1*Yb1*(-6.0*Yt*Yt+3.0*la3)+la1*la4+la2*la4+la4*la4+la5*la5))/(64.0*pi*pi);
106 //beta_lambda_4
107 beta[12] = (3.0*g1*g1*(g2*g2-la4)-9.0*g2*g2*la4+6.0*Yb2*Yb2*la4+6.0*Yt*Yt*la4+2.0*Ytau1*Ytau1*la4
108 +2.0*Ytau2*Ytau2*la4+2.0*la1*la4+2.0*la2*la4+8.0*la3*la4+4.0*la4*la4+6.0*Yb1*Yb1*(2.0*Yt*Yt+la4)
109 +8.0*la5*la5)/(16.0*pi*pi);
110 //beta_lambda_5
111 beta[13] = (-3.0*g1*g1-9.0*g2*g2
112 +2.0*(3.0*Yb1*Yb1+3.0*Yb2*Yb2+3.0*Yt*Yt+Ytau1*Ytau1+Ytau2*Ytau2+la1+la2+4.0*la3+6.0*la4))*la5/(16.0*pi*pi);
113
114 if(ord==1||ord==2){
115 //beta_g1
116 beta[0] += (g1*g1*g1*(88.0*g3*g3-17.0*Yt*Yt))/(1536.0*pi*pi*pi*pi);
117 //beta_g2
118 beta[1] += (3.0*g2*g2*g2*(8.0*g3*g3-Yt*Yt))/(512.0*pi*pi*pi*pi);
119 //beta_g3
120 beta[2] += -((g3*g3*g3*(13.0*g3*g3+Yt*Yt))/(128.0*pi*pi*pi*pi));
121 //beta_Yt
122 beta[3] += (Yt*(-216.0*g3*g3*g3*g3+72.0*g3*g3*Yt*Yt-24.0*Yt*Yt*Yt*Yt-12.0*Yt*Yt*la2
123 +3.0*la2*la2+2.0*la3*la3+2.0*la3*la4+2.0*la4*la4+3.0*la5*la5))/(512.0*pi*pi*pi*pi);
124 //beta_Yb
125 beta[4] += (-1296.0*g3*g3*g3*g3*(Yb1+Yb2)+16.0*g3*g3*(4.0*Yb1+3.0*Yb2)*Yt*Yt
126 -3.0*(2.0*Yb1*(5.0*Yt*Yt*Yt*Yt-3.0*la1*la1-2.0*la3*la3
127 +4.0*Yt*Yt*(la3-la4)-2.0*la3*la4
128 -2.0*la4*la4-3.0*la5*la5)
129 +Yb2*(Yt*Yt*Yt*Yt-2.0*(3.0*la2*la2+2.0*la3*la3+2.0*la3*la4
130 +2.0*la4*la4+3.0*la5*la5))))
131 /(3072.0*pi*pi*pi*pi);
132 //beta_Ytau
133 beta[5] += (80.0*g3*g3*Yt*Yt*Ytau2-27.0*Yt*Yt*Yt*Yt*Ytau2+6.0*Ytau1*la1*la1+6.0*Ytau2*la2*la2
134 +(Ytau1+Ytau2)*(4.0*la3*la3+4.0*la3*la4+4.0*la4*la4+6.0*la5*la5))
135 /(1024.0*pi*pi*pi*pi);
136 //beta_m11_2
137 beta[6] += -(8.0*m22_2*(2.0*la3*la3+2.0*la3*la4+2.0*la4*la4
138 +3.0*Yt*Yt*(2.0*la3+la4)+3.0*la5*la5)
139 +m11_2*(30.0*la1*la1+4.0*la3*la3+4.0*la3*la4
140 +4.0*la4*la4+6.0*la5*la5))/(512.0*pi*pi*pi*pi);
141 //beta_m22_2
142 beta[7] += -((-80.0*g3*g3*m22_2*Yt*Yt
143 +8.0*m11_2*(2.0*la3*la3+2.0*la3*la4+2.0*la4*la4+3.0*la5*la5)
144 +m22_2*(27.0*Yt*Yt*Yt*Yt+72.0*Yt*Yt*la2+30.0*la2*la2+4.0*la3*la3
145 +4.0*la3*la4+4.0*la4*la4+6.0*la5*la5))/(512.0*pi*pi*pi*pi));
146 //beta_m12_2
147 beta[8] += (m12_2*(72.0*(la1*la1+la2*la2-4.0*la3*la4-8.0*la3*la5-8.0*la4*la5
148 +2.0*la5*la5-4.0*(la1+la2)*(la3+la4+la5))
149 +Yt*Yt*(960.0*g3*g3
150 -36.0*(9.0*Yt*Yt+8.0*(la3+2.0*la4+3.0*la5)))))
151 /(12288.0*pi*pi*pi*pi);
152 //beta_lambda_1
153 beta[9] += -(39.0*la1*la1*la1
154 +la1*(10.0*la3*la3+10.0*la3*la4+6.0*la4*la4+7.0*la5*la5)
155 +2.0*(4.0*la3*la3*la3+6.0*la3*la3*la4+8.0*la3*la4*la4
156 +3.0*la4*la4*la4+10.0*la3*la5*la5+11.0*la4*la5*la5
157 +3.0*Yt*Yt*(2.0*la3*la3+2.0*la3*la4+la4*la4+la5*la5)))
158 /(128.0*pi*pi*pi*pi);
159 //beta_lambda_2
160 beta[10] += -(-60.0*pow(Yt,6)+3.0*Yt*Yt*Yt*Yt*la2+72.0*Yt*Yt*la2*la2
161 +16.0*g3*g3*(4.0*Yt*Yt*Yt*Yt-5.0*Yt*Yt*la2)
162 +2.0*(39.0*la2*la2*la2+10.0*la2*la3*la3+8.0*la3*la3*la3
163 +10.0*la2*la3*la4+12.0*la3*la3*la4+6.0*la2*la4*la4
164 +16.0*la3*la4*la4+6.0*la4*la4*la4+7.0*la2*la5*la5
165 +20.0*la3*la5*la5+22.0*la4*la5*la5))/(256.0*pi*pi*pi*pi);
166 //beta_lambda_3
167 beta[11] += -(-80.0*g3*g3*Yt*Yt*la3+27.0*Yt*Yt*Yt*Yt*la3
168 +12.0*Yt*Yt*(2.0*la3*la3+la4*la4+2.0*la2*(3.0*la3+la4)+la5*la5)
169 +2.0*(la1*la1*(15.0*la3+4.0*la4)+la2*la2*(15.0*la3+4.0*la4)
170 +2.0*(la1+la2)*(18.0*la3*la3+8.0*la3*la4+7.0*la4*la4+9.0*la5*la5)
171 +2.0*(6.0*la3*la3*la3+2.0*la3*la3*la4+8.0*la3*la4*la4
172 +6.0*la4*la4*la4+9.0*la3*la5*la5+22.0*la4*la5*la5)))
173 /(512.0*pi*pi*pi*pi);
174 //beta_lambda_4
175 beta[12] += -(-80.0*g3*g3*Yt*Yt*la4+27.0*Yt*Yt*Yt*Yt*la4
176 +24.0*Yt*Yt*(la2*la4+2.0*la3*la4+la4*la4+2.0*la5*la5)
177 +2.0*(7.0*la1*la1*la4+7.0*la2*la2*la4+28.0*la3*la3*la4
178 +28.0*la3*la4*la4+48.0*la3*la5*la5+26.0*la4*la5*la5
179 +4.0*(la1+la2)*(10.0*la3*la4+5.0*la4*la4+6.0*la5*la5))
180 )/(512.0*pi*pi*pi*pi);
181 //beta_lambda_5
182 beta[13] += (la5*(80.0*g3*g3*Yt*Yt-3.0*Yt*Yt*Yt*Yt-24.0*Yt*Yt*(la2+2.0*la3+3.0*la4)
183 -2.0*(7.0*la1*la1+7.0*la2*la2+40.0*la1*la3+40.0*la2*la3+28.0*la3*la3
184 +44.0*la1*la4+44.0*la2*la4+76.0*la3*la4+32.0*la4*la4-6.0*la5*la5))
185 )/(512 *pi*pi*pi*pi);
186
187 if(ord==2){
188 //beta_g1
189 beta[0] += (g1*g1*g1*(208.0*g1*g1+108.0*g2*g2-15.0*Yb1*Yb1-15.0*Yb2*Yb2
190 -45.0*Ytau1*Ytau1-45.0*Ytau2*Ytau2))/(4608.0*pi*pi*pi*pi);
191 //beta_g2
192 beta[1] += (g2*g2*g2*(4.0*g1*g1+16.0*g2*g2-3.0*Yb1*Yb1-3.0*Yb2*Yb2
193 -Ytau1*Ytau1-Ytau2*Ytau2))/(512.0*pi*pi*pi*pi);
194 //beta_g3
195 beta[2] += (g3*g3*g3*(11.0*g1*g1+3.0*(9.0*g2*g2-4.0*(Yb1*Yb1+Yb2*Yb2))))/(1536.0*pi*pi*pi*pi);
196 //beta_Yt
197 beta[3] += (2534.0*g1*g1*g1*g1
198 +g1*g1*(-324.0*g2*g2+912.0*g3*g3-123.0*Yb1*Yb1+63.0*Yb2*Yb2
199 +3537.0*Yt*Yt+1350.0*Ytau2*Ytau2)
200 -9.0*(252.0*g2*g2*g2*g2
201 -9.0*g2*g2*(48.0*g3*g3+11.0*Yb1*Yb1+33.0*Yb2*Yb2+75.0*Yt*Yt+10.0*Ytau2*Ytau2)
202 -4.0*(16.0*g3*g3*(4.0*Yb1*Yb1+3.0*Yb2*Yb2)
203 -3.0*(10.0*Yb1*Yb1*Yb1*Yb1+Yb2*Yb2*Yb2*Yb2
204 +Yb2*Yb2*(11.0*Yt*Yt-5.0*Ytau2*Ytau2)
205 +9.0*Ytau2*Ytau2*(Yt*Yt+Ytau2*Ytau2)
206 +Yb1*Yb1*(10.0*Yt*Yt+3.0*Ytau1*Ytau1+8.0*la3-8.0*la4)))))
207 *Yt/(110592.0*pi*pi*pi*pi);
208 //beta_Yb
209 beta[4] += -(226.0*g1*g1*g1*g1*(Yb1+Yb2)
210 +3.0*g1*g1*(-711.0*Yb1*Yb1*Yb1-711.0*Yb2*Yb2*Yb2+324.0*g2*g2*(Yb1+Yb2)
211 -496.0*g3*g3*(Yb1+Yb2)+53.0*Yb1*Yt*Yt-273.0*Yb2*Yt*Yt
212 -450.0*Yb1*Ytau1*Ytau1-450.0*Yb2*Ytau2*Ytau2)
213 +27.0*(84.0*g2*g2*g2*g2*(Yb1+Yb2)
214 -3.0*g2*g2*(75.0*Yb1*Yb1*Yb1+75.0*Yb2*Yb2*Yb2+48.0*g3*g3*(Yb1+Yb2)
215 +11.0*Yb1*Yt*Yt+33.0*Yb2*Yt*Yt+10.0*Yb1*Ytau1*Ytau1
216 +10.0*Yb2*Ytau2*Ytau2)
217 +4.0*(48.0*pow(Yb1,5)-144.0*g3*g3*(Yb1*Yb1*Yb1+Yb2*Yb2*Yb2)
218 +3.0*Yb1*(3.0*Ytau1*Ytau1*Ytau1*Ytau1+Yt*Yt*Ytau2*Ytau2)
219 +Yb1*Yb1*Yb1*(10.0*Yt*Yt+9.0*Ytau1*Ytau1+24.0*la1)
220 +Yb2*(48.0*Yb2*Yb2*Yb2*Yb2-5.0*Yt*Yt*Ytau2*Ytau2
221 +9.0*Ytau2*Ytau2*Ytau2*Ytau2
222 +Yb2*Yb2*(11.0*Yt*Yt+9.0*Ytau2*Ytau2+24.0*la2)))))
223 /(110592.0*pi*pi*pi*pi);
224 //beta_Ytau
225 beta[5] += (966.0*g1*g1*g1*g1*(Ytau1+Ytau2)
226 +g1*g1*(50.0*Yb1*Yb1*Ytau1+537.0*Ytau1*Ytau1*Ytau1+50.0*Yb2*Yb2*Ytau2
227 +170.0*Yt*Yt*Ytau2+537.0*Ytau2*Ytau2*Ytau2+108.0*g2*g2*(Ytau1+Ytau2))
228 -3.0*(84.0*g2*g2*g2*g2*(Ytau1+Ytau2)
229 -15.0*g2*g2*(6.0*Yb1*Yb1*Ytau1+11.0*Ytau1*Ytau1*Ytau1+6.0*Yb2*Yb2*Ytau2
230 +6.0*Yt*Yt*Ytau2+11.0*Ytau2*Ytau2*Ytau2)
231 +4.0*(-80.0*g3*g3*(Yb1*Yb1*Ytau1+Yb2*Yb2*Ytau2)
232 +3.0*(9.0*Yb1*Yb1*Yb1*Yb1*Ytau1+4.0*pow(Ytau1,5)+9.0*Yb2*Yb2*Yb2*Yb2*Ytau2
233 -2.0*Yb2*Yb2*Yt*Yt*Ytau2+9.0*Yb2*Yb2*Ytau2*Ytau2*Ytau2
234 +9.0*Yt*Yt*Ytau2*Ytau2*Ytau2+4.0*pow(Ytau2,5)
235 +3.0*Yb1*Yb1*(3.0*Ytau1*Ytau1*Ytau1+Yt*Yt*(Ytau1+Ytau2))
236 +8.0*Ytau1*Ytau1*Ytau1*la1+8.0*Ytau2*Ytau2*Ytau2*la2))))
237 /(12288.0*pi*pi*pi*pi);
238 //beta_m11_2
239 beta[6] += (3.0*g1*g1*g1*g1*(193.0*m11_2+40.0*m22_2)
240 +2.0*g1*g1*(45.0*g2*g2*m11_2
241 +2.0*(m11_2*(25.0*Yb1*Yb1+75.0*Ytau1*Ytau1+144.0*la1)
242 +48.0*m22_2*(2.0*la3+la4)))
243 -3.0*(3.0*g2*g2*g2*g2*(41.0*m11_2-40.0*m22_2)
244 -12.0*g2*g2*(m11_2*(15.0*Yb1*Yb1+5.0*Ytau1*Ytau1+48.0*la1)
245 +16.0*m22_2*(2.0*la3+la4))
246 +8.0*(-80.0*g3*g3*m11_2*Yb1*Yb1
247 +3.0*m11_2*(9.0*Yb1*Yb1*Yb1*Yb1+3.0*Ytau1*Ytau1*Ytau1*Ytau1
248 +8.0*Ytau1*Ytau1*la1+3.0*Yb1*Yb1*(Yt*Yt+8.0*la1))
249 +8.0*m22_2*(3.0*Yb2*Yb2+Ytau2*Ytau2)*(2.0*la3+la4))))
250 /(12288.0*pi*pi*pi*pi);
251 //beta_m22_2
252 beta[7] += (3.0*g1*g1*g1*g1*(40.0*m11_2+193.0*m22_2)
253 +2.0*g1*g1*(45.0*g2*g2*m22_2
254 +2.0*(m22_2*(25.0*Yb2*Yb2+85.0*Yt*Yt+75.0*Ytau2*Ytau2+144.0*la2)
255 +48.0*m11_2*(2.0*la3+la4)))
256 +3.0*(3.0*g2*g2*g2*g2*(40.0*m11_2-41.0*m22_2)
257 +12.0*g2*g2*(m22_2*(15.0*Yb2*Yb2+15.0*Yt*Yt+5.0*Ytau2*Ytau2+48.0*la2)
258 +16.0*m11_2*(2.0*la3+la4))
259 +8.0*(80.0*g3*g3*m22_2*Yb2*Yb2
260 -3.0*m22_2*(9.0*Yb2*Yb2*Yb2*Yb2+3.0*Yb1*Yb1*Yt*Yt
261 +3.0*Ytau2*Ytau2*Ytau2*Ytau2+8.0*Ytau2*Ytau2*la2
262 +2.0*Yb2*Yb2*(7.0*Yt*Yt+12.0*la2))
263 -8.0*m11_2*(3.0*Yb1*Yb1+Ytau1*Ytau1)*(2.0*la3+la4))))
264 /(12288.0*pi*pi*pi*pi);
265 //beta_m12_2
266 beta[8] += (9.0*(51.0*g1*g1*g1*g1+10.0*g1*g1*g2*g2-81.0*g2*g2*g2*g2)
267 -324.0*(Yb1*Yb1*Yb1*Yb1+Yb2*Yb2*Yb2*Yb2)
268 +2.0*(85.0*g1*g1+135.0*g2*g2-432.0*Yb1*Yb1)*Yt*Yt
269 -108.0*(Ytau1*Ytau1*Ytau1*Ytau1+Ytau2*Ytau2*Ytau2*Ytau2)
270 +2.0*(Yb1*Yb1+Yb2*Yb2)*(25.0*g1*g1
271 +3.0*(45.0*g2*g2+4.0*(40.0*g3*g3+3.0*Yt*Yt-12.0*la3
272 -24.0*la4-36.0*la5)))
273 +192.0*(g1*g1+3.0*g2*g2)*(la3+2.0*la4+3.0*la5)
274 +6.0*(Ytau1*Ytau1+Ytau2*Ytau2)*(25.0*g1*g1+15.0*g2*g2
275 -16.0*(la3+2.0*la4+3.0*la5)))
276 *m12_2/(12288.0*pi*pi*pi*pi);
277 //beta_lambda_1
278 beta[9] += (-393.0*pow(g1,6)
279 +g1*g1*g1*g1*(-573.0*g2*g2+60.0*Yb1*Yb1-300.0*Ytau1*Ytau1
280 +651.0*la1+120.0*la3+60.0*la4)
281 +3.0*(291.0*pow(g2,6)
282 -3.0*g2*g2*g2*g2*(12.0*Yb1*Yb1+4.0*Ytau1*Ytau1+17.0*la1-40.0*la3-20.0*la4)
283 +12.0*g2*g2*(15.0*Yb1*Yb1*la1+5.0*Ytau1*Ytau1*la1
284 +4.0*(9.0*la1*la1+(2.0*la3+la4)*(2.0*la3+la4)))
285 -8.0*(-60.0*pow(Yb1,6)-20.0*pow(Ytau1,6)+Ytau1*Ytau1*Ytau1*Ytau1*la1
286 +24.0*Ytau1*Ytau1*la1*la1+3.0*Yb1*Yb1*Yb1*Yb1*(-4.0*Yt*Yt+la1)
287 +9.0*Yb1*Yb1*la1*(Yt*Yt+8.0*la1)
288 +16.0*g3*g3*(4.0*Yb1*Yb1*Yb1*Yb1-5.0*Yb1*Yb1*la1)
289 +24.0*Yb2*Yb2*la3*la3+8.0*Ytau2*Ytau2*la3*la3+24.0*Yb2*Yb2*la3*la4
290 +8.0*Ytau2*Ytau2*la3*la4+12.0*Yb2*Yb2*la4*la4+4.0*Ytau2*Ytau2*la4*la4
291 +12.0*Yb2*Yb2*la5*la5+4.0*Ytau2*Ytau2*la5*la5))
292 +g1*g1*(-303.0*g2*g2*g2*g2
293 +6.0*g2*g2*(36.0*Yb1*Yb1+44.0*Ytau1*Ytau1+39.0*la1+20.0*la4)
294 +4.0*(16.0*Yb1*Yb1*Yb1*Yb1+25.0*Yb1*Yb1*la1
295 +3.0*(-16.0*Ytau1*Ytau1*Ytau1*Ytau1+25.0*Ytau1*Ytau1*la1+36.0*la1*la1
296 +16.0*la3*la3+16.0*la3*la4+8.0*la4*la4-4.0*la5*la5))))
297 /(6144.0*pi*pi*pi*pi);
298 //beta_lambda_2
299 beta[10] += (-393.0*pow(g1,6)
300 +g1*g1*g1*g1*(-573.0*g2*g2+60.0*Yb2*Yb2-228.0*Yt*Yt-300.0*Ytau2*Ytau2
301 +651.0*la2+120.0*la3+60.0*la4)
302 +g1*g1*(-303.0*g2*g2*g2*g2
303 +6.0*g2*g2*(36.0*Yb2*Yb2+84.0*Yt*Yt+44.0*Ytau2*Ytau2+39.0*la2+20.0*la4)
304 +4.0*(16.0*Yb2*Yb2*Yb2*Yb2-32.0*Yt*Yt*Yt*Yt-48.0*Ytau2*Ytau2*Ytau2*Ytau2
305 +25.0*Yb2*Yb2*la2+85.0*Yt*Yt*la2+75.0*Ytau2*Ytau2*la2
306 +108.0*la2*la2+48.0*la3*la3+48.0*la3*la4+24.0*la4*la4-12.0*la5*la5))
307 +3.0*(291.0*pow(g2,6)
308 -3.0*g2*g2*g2*g2*(12.0*Yb2*Yb2+12.0*Yt*Yt+4.0*Ytau2*Ytau2
309 +17.0*la2-40.0*la3-20.0*la4)
310 +12.0*g2*g2*(15.0*Yb2*Yb2*la2+15.0*Yt*Yt*la2+5.0*Ytau2*Ytau2*la2
311 +36.0*la2*la2+16.0*la3*la3+16.0*la3*la4+4.0*la4*la4)
312 -8.0*(-60.0*pow(Yb2,6)-12.0*Yb1*Yb1*Yt*Yt*Yt*Yt-20.0*pow(Ytau2,6)
313 +9.0*Yb1*Yb1*Yt*Yt*la2+Ytau2*Ytau2*Ytau2*Ytau2*la2
314 +24.0*Ytau2*Ytau2*la2*la2+3.0*Yb2*Yb2*Yb2*Yb2*(4.0*Yt*Yt+la2)
315 +16.0*g3*g3*(4.0*Yb2*Yb2*Yb2*Yb2-5.0*Yb2*Yb2*la2)
316 +6.0*Yb2*Yb2*(2.0*Yt*Yt*Yt*Yt+7.0*Yt*Yt*la2+12.0*la2*la2)
317 +24.0*Yb1*Yb1*la3*la3+8.0*Ytau1*Ytau1*la3*la3+24.0*Yb1*Yb1*la3*la4
318 +8.0*Ytau1*Ytau1*la3*la4+12.0*Yb1*Yb1*la4*la4+4.0*Ytau1*Ytau1*la4*la4
319 +12.0*Yb1*Yb1*la5*la5+4.0*Ytau1*Ytau1*la5*la5)))/(6144.0*pi*pi*pi*pi);
320 //beta_lambda_3
321 beta[11] += (-393.0*pow(g1,6)
322 +3.0*g1*g1*g1*g1*(101.0*g2*g2+10.0*Yb1*Yb1+10.0*Yb2*Yb2-38.0*Yt*Yt-50.0*Ytau1*Ytau1
323 -50.0*Ytau2*Ytau2+30.0*la1+30.0*la2+197.0*la3+20.0*la4)
324 +g1*g1*(33.0*g2*g2*g2*g2-6.0*g2*g2*(18.0*Yb1*Yb1+18.0*Yb2*Yb2+42.0*Yt*Yt
325 +22.0*Ytau1*Ytau1+22.0*Ytau2*Ytau2+10.0*la1
326 +10.0*la2-11.0*la3+12.0*la4)
327 +2.0*(25.0*Yb2*Yb2*la3+85.0*Yt*Yt*la3+75.0*Ytau1*Ytau1*la3
328 +75.0*Ytau2*Ytau2*la3+144.0*la1*la3+144.0*la2*la3+24*la3*la3
329 +Yb1*Yb1*(-16.0*Yt*Yt+25.0*la3)
330 +48.0*la1*la4+48.0*la2*la4-24.0*la4*la4+48.0*la5*la5))
331 +3.0*(291.0*pow(g2,6)-3.0*g2*g2*g2*g2*(6.0*Yb1*Yb1+6.0*Yb2*Yb2+6.0*Yt*Yt
332 +2.0*Ytau1*Ytau1+2.0*Ytau2*Ytau2
333 -30.0*la1-30.0*la2+37.0*la3-20.0*la4)
334 +6.0*g2*g2*(15.0*Yb1*Yb1*la3+15.0*Yb2*Yb2*la3+15.0*Yt*Yt*la3
335 +5.0*Ytau1*Ytau1*la3+5.0*Ytau2*Ytau2*la3+48.0*la1*la3
336 +48.0*la2*la3+8.0*la3*la3+24.0*la1*la4+24.0*la2*la4
337 -16.0*la3*la4+8.0*la4*la4)
338 -4.0*(-9.0*Yb1*Yb1*Yb1*Yb1*(8.0*Yt*Yt-3.0*la3)+27.0*Yb2*Yb2*Yb2*Yb2*la3
339 +42.0*Yb2*Yb2*Yt*Yt*la3+9.0*Ytau1*Ytau1*Ytau1*Ytau1*la3
340 +9.0*Ytau2*Ytau2*Ytau2*Ytau2*la3+24.0*Ytau1*Ytau1*la1*la3
341 +72.0*Yb2*Yb2*la2*la3+24.0*Ytau2*Ytau2*la2*la3+24.0*Yb2*Yb2*la3*la3
342 +8.0*Ytau1*Ytau1*la3*la3+8.0*Ytau2*Ytau2*la3*la3
343 +16.0*g3*g3*(Yb1*Yb1*(8.0*Yt*Yt-5.0*la3)-5.0*Yb2*Yb2*la3)
344 +48.0*Yb2*Yb2*Yt*Yt*la4+8.0*Ytau1*Ytau1*la1*la4+24.0*Yb2*Yb2*la2*la4
345 +8.0*Ytau2*Ytau2*la2*la4+12.0*Yb2*Yb2*la4*la4+4.0*Ytau1*Ytau1*la4*la4
346 +4.0*Ytau2*Ytau2*la4*la4+12.0*Yb2*Yb2*la5*la5+4.0*Ytau1*Ytau1*la5*la5
347 +4.0*Ytau2*Ytau2*la5*la5
348 -6.0*Yb1*Yb1*(12.0*Yt*Yt*Yt*Yt+5.0*Yt*Yt*la3
349 -2.0*(6.0*la1*la3+2.0*la3*la3+2.0*la1*la4+la4*la4+la5*la5)))))
350 /(6144.0*pi*pi*pi*pi);
351 //beta_lambda_4
352 beta[12] += (g1*g1*g1*g1*(-876.0*g2*g2+471.0*la4)
353 +2.0*g1*g1*(-168.0*g2*g2*g2*g2+25.0*Yb2*Yb2*la4+85.0*Yt*Yt*la4
354 +75.0*Ytau1*Ytau1*la4+75.0*Ytau2*Ytau2*la4
355 +48.0*la1*la4+48.0*la2*la4+48.0*la3*la4+96.0*la4*la4
356 +Yb1*Yb1*(16.0*Yt*Yt+25.0*la4)
357 +3.0*g2*g2*(36.0*Yb1*Yb1+36.0*Yb2*Yb2+84.0*Yt*Yt
358 +44.0*Ytau1*Ytau1+44.0*Ytau2*Ytau2
359 +20.0*la1+20.0*la2+8.0*la3+51.0*la4)
360 +192.0*la5*la5)
361 -3.0*(231.0*g2*g2*g2*g2*la4-90.0*g2*g2*Yb2*Yb2*la4
362 +108.0*Yb2*Yb2*Yb2*Yb2*la4-90.0*g2*g2*Yt*Yt*la4
363 -216.0*Yb2*Yb2*Yt*Yt*la4-30.0*g2*g2*Ytau1*Ytau1*la4
364 +36.0*Ytau1*Ytau1*Ytau1*Ytau1*la4-30.0*g2*g2*Ytau2*Ytau2*la4
365 +36.0*Ytau2*Ytau2*Ytau2*Ytau2*la4+32.0*Ytau1*Ytau1*la1*la4
366 +96.0*Yb2*Yb2*la2*la4+32.0*Ytau2*Ytau2*la2*la4-288.0*g2*g2*la3*la4
367 +192.0*Yb2*Yb2*la3*la4+64.0*Ytau1*Ytau1*la3*la4+64.0*Ytau2*Ytau2*la3*la4
368 -144.0*g2*g2*la4*la4+96.0*Yb2*Yb2*la4*la4+32.0*Ytau1*Ytau1*la4*la4
369 +32.0*Ytau2*Ytau2*la4*la4+12.0*Yb1*Yb1*Yb1*Yb1*(16.0*Yt*Yt+9.0*la4)
370 -64.0*g3*g3*(5.0*Yb2*Yb2*la4+Yb1*Yb1*(8.0*Yt*Yt+5.0*la4))
371 -432.0*g2*g2*la5*la5+192.0*Yb2*Yb2*la5*la5+64.0*Ytau1*Ytau1*la5*la5
372 +64.0*Ytau2*Ytau2*la5*la5
373 +6.0*Yb1*Yb1*(32.0*Yt*Yt*Yt*Yt-15.0*g2*g2*la4
374 +4.0*Yt*Yt*(8.0*la3+11.0*la4)
375 +16.0*(la1*la4+2.0*la3*la4+la4*la4+2.0*la5*la5))))
376 /(6144.0*pi*pi*pi*pi);
377 //beta_lambda_5
378 beta[13] += (471.0*g1*g1*g1*g1
379 +2.0*g1*g1*(57.0*g2*g2+25.0*Yb1*Yb1+25.0*Yb2*Yb2+85.0*Yt*Yt+75.0*Ytau1*Ytau1
380 +75.0*Ytau2*Ytau2-24.0*la1-24.0*la2+192.0*la3+288.0*la4)
381 -3.0*(231.0*g2*g2*g2*g2
382 -6.0*g2*g2*(15.0*Yb1*Yb1+15.0*Yb2*Yb2+15.0*Yt*Yt+5.0*Ytau1*Ytau1
383 +5.0*Ytau2*Ytau2+48.0*la3+96.0*la4)
384 +4.0*(3.0*Yb1*Yb1*Yb1*Yb1+3.0*Yb2*Yb2*Yb2*Yb2-80.0*g3*g3*(Yb1*Yb1+Yb2*Yb2)
385 -6.0*Yb2*Yb2*Yt*Yt+Ytau1*Ytau1*Ytau1*Ytau1+Ytau2*Ytau2*Ytau2*Ytau2
386 +8.0*Ytau1*Ytau1*la1+24.0*Yb2*Yb2*la2+8.0*Ytau2*Ytau2*la2+48.0*Yb2*Yb2*la3
387 +16.0*Ytau1*Ytau1*la3+16.0*Ytau2*Ytau2*la3+72.0*Yb2*Yb2*la4
388 +24.0*Ytau1*Ytau1*la4+24.0*Ytau2*Ytau2*la4
389 +6.0*Yb1*Yb1*(11.0*Yt*Yt+4.0*la1+8.0*la3+12.0*la4))))
390 *la5/(6144.0*pi*pi*pi*pi);
391 }
392 }
393
394 return 0;
395}
396
397int JacobianZ2 (double t, const double y[], double *dfdy, double dfdt[], void *order)
398{
399 return 0;
400}
401
402int RGEcheckZ2(const double InitialValues[], const double t1, const double Rpeps, const double tNLOuni)
403{
404 int check=0;
405
406 //perturbativity of the Yukawa couplings
407 for(int i=3;i<6;i++)
408 {
409 if(fabs(InitialValues[i])>sqrt(4.0*M_PI)) check=1;
410 }
411
412 //perturbativity of the quartic Higgs couplings
413 for(int i=9;i<14;i++)
414 {
415 if(fabs(InitialValues[i])>4.0*M_PI) check=1;
416 }
417
418 //positivity
419 if(InitialValues[9]<0.0) check=1;
420 if(InitialValues[10]<0.0) check=1;
421 if(InitialValues[11]+sqrt(fabs(InitialValues[9]*InitialValues[10]))<0.0) check=1;
422 if(InitialValues[11]+InitialValues[12]-fabs(InitialValues[13])+sqrt(fabs(InitialValues[9]*InitialValues[10]))<0.0) check=1;
423
424 //NLO unitarity
425 double YtQ = InitialValues[3];
426 if(t1>tNLOuni)
427 {
428 double la1Q = InitialValues[9];
429 double la2Q = InitialValues[10];
430 double la3Q = InitialValues[11];
431 double la4Q = InitialValues[12];
432 double la5Q = InitialValues[13];
433
434 double betalambda1 = 6.0*la1Q*la1Q + 2.0*la3Q*la3Q + 2.0*la3Q*la4Q + la4Q*la4Q + la5Q*la5Q;
435// + 6.0*la1Q*Yb1Q*Yb1Q + 2.0*la1Q*Ytau1Q*Ytau1Q
436// - 6.0*Yb1Q*Yb1Q*Yb1Q*Yb1Q - 2.0*Ytau1Q*Ytau1Q*Ytau1Q*Ytau1Q;
437 double betalambda2 = 6.0*la2Q*la2Q + 2.0*la3Q*la3Q + 2.0*la3Q*la4Q + la4Q*la4Q + la5Q*la5Q + 6.0*la2Q*YtQ*YtQ - 6.0*YtQ*YtQ*YtQ*YtQ;
438// + 6.0*la2Q*Yb2Q*Yb2Q + 2.0*la2Q*Ytau2Q*Ytau2Q
439// - 6.0*Yb2Q*Yb2Q*Yb2Q*Yb2Q - 2.0*Ytau2Q*Ytau2Q*Ytau2Q*Ytau2Q;
440 double betalambda3 = 3.0*la1Q*la3Q + 3.0*la2Q*la3Q + 2.0*la3Q*la3Q + la1Q*la4Q + la2Q*la4Q + la4Q*la4Q + la5Q*la5Q + 3.0*la3Q*YtQ*YtQ;
441// + 3.0*la3Q*Yb1Q*Yb1Q + 3.0*la3Q*Yb2Q*Yb2Q + la3Q*Ytau1Q*Ytau1Q + la3Q*Ytau2Q*Ytau2Q
442// - 6.0*Yb1Q*Yb1Q*(Yb2Q*Yb2Q + YtQ*YtQ) - 2.0*Ytau1Q*Ytau1Q*Ytau2Q*Ytau2Q;
443 double betalambda4 = la1Q*la4Q + la2Q*la4Q + 4.0*la3Q*la4Q + 2.0*la4Q*la4Q + 4.0*la5Q*la5Q + 3.0*la4Q*YtQ*YtQ;
444// + 3.0*la4Q*Yb1Q*Yb1Q + 3.0*la4Q*Yb2Q*Yb2Q + la4Q*Ytau1Q*Ytau1Q + la4Q*Ytau2Q*Ytau2Q
445// - 2.0*Ytau1Q*Ytau1Q*Ytau2Q*Ytau2Q - 6.0*Yb1Q*Yb1Q*(Yb2Q*Yb2Q - YtQ*YtQ);
446 double betalambda5 = la5Q*(la1Q + la2Q + 4.0*la3Q + 6.0*la4Q + 3.0*YtQ*YtQ);
447// + 3.0*la5Q*Yb1Q*Yb1Q + la5Q*Ytau1Q*Ytau1Q + la5Q*(3.0*Yb2Q*Yb2Q + Ytau2Q*Ytau2Q)
448// - 6.0*Yb1Q*Yb1Q*Yb2Q*Yb2Q - 2.0*Ytau1Q*Ytau1Q*Ytau2Q*Ytau2Q;
449
450 double uniLO1 = -3.0*la1Q/(16.0*M_PI);
451 gslpp::complex uniNLO1 = -3.0*la1Q/(16.0*M_PI) +9.0*betalambda1/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(9.0*la1Q*la1Q+(2.0*la3Q+la4Q)*(2.0*la3Q+la4Q))/(256.0*M_PI*M_PI*M_PI);
452 double uniLO2 = -3.0*la2Q/(16.0*M_PI);
453 gslpp::complex uniNLO2 = -3.0*la2Q/(16.0*M_PI) +9.0*betalambda2/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(9.0*la2Q*la2Q+(2.0*la3Q+la4Q)*(2.0*la3Q+la4Q))/(256.0*M_PI*M_PI*M_PI);
454 double uniLO3 = -(2.0*la3Q+la4Q)/(16.0*M_PI);
455 gslpp::complex uniNLO3 = -(2.0*la3Q+la4Q)/(16.0*M_PI) +3.0*(2.0*betalambda3+betalambda4)/(256.0*M_PI*M_PI*M_PI) +3.0*(gslpp::complex::i()*M_PI-1.0)*(la1Q+la2Q)*(2.0*la3Q+la4Q)/(256.0*M_PI*M_PI*M_PI);
456 double uniLO4 = -(la3Q+2.0*la4Q)/(16.0*M_PI);
457 gslpp::complex uniNLO4 = -(la3Q+2.0*la4Q)/(16.0*M_PI) +3.0*(betalambda3+2.0*betalambda4)/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la3Q*la3Q+4.0*la3Q*la4Q+4.0*la4Q*la4Q+9.0*la5Q*la5Q)/(256.0*M_PI*M_PI*M_PI);
458 double uniLO5 = -3.0*la5Q/(16.0*M_PI);
459 gslpp::complex uniNLO5 = -3.0*la5Q/(16.0*M_PI) +9.0*betalambda5/(256.0*M_PI*M_PI*M_PI) +3.0*(gslpp::complex::i()*M_PI-1.0)*(la3Q+2.0*la4Q)*la5Q/(128.0*M_PI*M_PI*M_PI);
460 double uniLO6 = -la1Q/(16.0*M_PI);
461 gslpp::complex uniNLO6 = -la1Q/(16.0*M_PI) +3.0*betalambda1/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la1Q*la1Q+la4Q*la4Q)/(256.0*M_PI*M_PI*M_PI);
462 double uniLO7 = -la2Q/(16.0*M_PI);
463 gslpp::complex uniNLO7 = -la2Q/(16.0*M_PI) +3.0*betalambda2/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la2Q*la2Q+la4Q*la4Q)/(256.0*M_PI*M_PI*M_PI);
464 double uniLO8 = -la4Q/(16.0*M_PI);
465 gslpp::complex uniNLO8 = -la4Q/(16.0*M_PI) +3.0*betalambda4/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la1Q+la2Q)*la4Q/(256.0*M_PI*M_PI*M_PI);
466 double uniLO10 = -la3Q/(16.0*M_PI);
467 gslpp::complex uniNLO10 = -la3Q/(16.0*M_PI) +3.0*betalambda3/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la3Q*la3Q+la5Q*la5Q)/(256.0*M_PI*M_PI*M_PI);
468 double uniLO11 = -la5Q/(16.0*M_PI);
469 gslpp::complex uniNLO11 = -la5Q/(16.0*M_PI) +3.0*betalambda5/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*la3Q*la5Q/(128.0*M_PI*M_PI*M_PI);
470 double uniLO14 = -(la3Q-la4Q)/(16.0*M_PI);
471 gslpp::complex uniNLO14 = -(la3Q-la4Q)/(16.0*M_PI) +3.0*(betalambda3-betalambda4)/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la3Q-la4Q)*(la3Q-la4Q)/(256.0*M_PI*M_PI*M_PI);
472 double uniLO15 = -la1Q/(16.0*M_PI);
473 gslpp::complex uniNLO15 = -la1Q/(16.0*M_PI) +3.0*betalambda1/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la1Q*la1Q+la5Q*la5Q)/(256.0*M_PI*M_PI*M_PI);
474 double uniLO16 = -la2Q/(16.0*M_PI);
475 gslpp::complex uniNLO16 = -la2Q/(16.0*M_PI) +3.0*betalambda2/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la2Q*la2Q+la5Q*la5Q)/(256.0*M_PI*M_PI*M_PI);
476 double uniLO17 = -la5Q/(16.0*M_PI);
477 gslpp::complex uniNLO17 = -la5Q/(16.0*M_PI) +3.0*betalambda5/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la1Q+la2Q)*la5Q/(256.0*M_PI*M_PI*M_PI);
478 double uniLO24 = -(la3Q+la4Q)/(16.0*M_PI);
479 gslpp::complex uniNLO24 = -(la3Q+la4Q)/(16.0*M_PI) +3.0*(betalambda3+betalambda4)/(256.0*M_PI*M_PI*M_PI) +(gslpp::complex::i()*M_PI-1.0)*(la3Q+la4Q)*(la3Q+la4Q)/(256.0*M_PI*M_PI*M_PI);
480
481 double uniLOev1 = 0.5*(uniLO1+uniLO2+sqrt((uniLO1-uniLO2)*(uniLO1-uniLO2)+4.0*uniLO3*uniLO3));
482 double uniLOev2 = 0.5*(uniLO1+uniLO2-sqrt((uniLO1-uniLO2)*(uniLO1-uniLO2)+4.0*uniLO3*uniLO3));
483 double uniLOev3 = uniLO4+uniLO5;
484 double uniLOev4 = uniLO4-uniLO5;
485 double uniLOev5 = 0.5*(uniLO6+uniLO7+sqrt((uniLO6-uniLO7)*(uniLO6-uniLO7)+4.0*uniLO8*uniLO8));
486 double uniLOev6 = 0.5*(uniLO6+uniLO7-sqrt((uniLO6-uniLO7)*(uniLO6-uniLO7)+4.0*uniLO8*uniLO8));
487 double uniLOev9 = uniLO10+uniLO11;
488 double uniLOev10 = uniLO10-uniLO11;
489 double uniLOev13 = 0.5*(uniLO15+uniLO16+sqrt((uniLO15-uniLO16)*(uniLO15-uniLO16)+4.0*uniLO17*uniLO17));
490 double uniLOev14 = 0.5*(uniLO15+uniLO16-sqrt((uniLO15-uniLO16)*(uniLO15-uniLO16)+4.0*uniLO17*uniLO17));
491 gslpp::complex uniNLOev1wowfr = 0.5*(uniNLO1+uniNLO2+sqrt((uniNLO1-uniNLO2)*(uniNLO1-uniNLO2)+4.0*uniNLO3*uniNLO3));
492 gslpp::complex uniNLOev2wowfr = 0.5*(uniNLO1+uniNLO2-sqrt((uniNLO1-uniNLO2)*(uniNLO1-uniNLO2)+4.0*uniNLO3*uniNLO3));
493 gslpp::complex uniNLOev3wowfr = uniNLO4+uniNLO5;
494 gslpp::complex uniNLOev4wowfr = uniNLO4-uniNLO5;
495 gslpp::complex uniNLOev5wowfr = 0.5*(uniNLO6+uniNLO7+sqrt((uniNLO6-uniNLO7)*(uniNLO6-uniNLO7)+4.0*uniNLO8*uniNLO8));
496 gslpp::complex uniNLOev6wowfr = 0.5*(uniNLO6+uniNLO7-sqrt((uniNLO6-uniNLO7)*(uniNLO6-uniNLO7)+4.0*uniNLO8*uniNLO8));
497 gslpp::complex uniNLOev9wowfr = uniNLO10+uniNLO11;
498 gslpp::complex uniNLOev10wowfr = uniNLO10-uniNLO11;
499 gslpp::complex uniNLOev13wowfr = 0.5*(uniNLO15+uniNLO16+sqrt((uniNLO15-uniNLO16)*(uniNLO15-uniNLO16)+4.0*uniNLO17*uniNLO17));
500 gslpp::complex uniNLOev14wowfr = 0.5*(uniNLO15+uniNLO16-sqrt((uniNLO15-uniNLO16)*(uniNLO15-uniNLO16)+4.0*uniNLO17*uniNLO17));
501
502 if( (uniNLOev1wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
503 if( (uniNLOev2wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
504 if( (uniNLOev3wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
505 if( (uniNLOev4wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
506 if( (uniNLOev5wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
507 if( (uniNLOev6wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
508 if( (uniNLOev9wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
509 if( (uniNLOev10wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
510 if( (uniNLOev13wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
511 if( (uniNLOev14wowfr - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
512 if( (uniNLO14 - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
513 if( (uniNLO24 - 0.5*gslpp::complex::i()).abs() > 0.5) check=1;
514
515 if( fabs(uniLOev1) > Rpeps && (uniNLOev1wowfr/uniLOev1-1.0).abs() > 1.0) check=1;
516 if( fabs(uniLOev2) > Rpeps && (uniNLOev2wowfr/uniLOev2-1.0).abs() > 1.0) check=1;
517 if( fabs(uniLOev3) > Rpeps && (uniNLOev3wowfr/uniLOev3-1.0).abs() > 1.0) check=1;
518 if( fabs(uniLOev4) > Rpeps && (uniNLOev4wowfr/uniLOev4-1.0).abs() > 1.0) check=1;
519 if( fabs(uniLOev5) > Rpeps && (uniNLOev5wowfr/uniLOev5-1.0).abs() > 1.0) check=1;
520 if( fabs(uniLOev6) > Rpeps && (uniNLOev6wowfr/uniLOev6-1.0).abs() > 1.0) check=1;
521 if( fabs(uniLOev9) > Rpeps && (uniNLOev9wowfr/uniLOev9-1.0).abs() > 1.0) check=1;
522 if( fabs(uniLOev10) > Rpeps && (uniNLOev10wowfr/uniLOev10-1.0).abs() > 1.0) check=1;
523 if( fabs(uniLOev13) > Rpeps && (uniNLOev13wowfr/uniLOev13-1.0).abs() > 1.0) check=1;
524 if( fabs(uniLOev14) > Rpeps && (uniNLOev14wowfr/uniLOev14-1.0).abs() > 1.0) check=1;
525 if( fabs(uniLO14) > Rpeps && (uniNLO14/uniLO14-1.0).abs() > 1.0) check=1;
526 if( fabs(uniLO24) > Rpeps && (uniNLO24/uniLO24-1.0).abs() > 1.0) check=1;
527 }
528
529 return check;
530}
531
532double GeneralTHDMZ2Runner::RGEGeneralTHDMZ2Runner(double InitialValues[], unsigned long int NumberOfRGEs, double Q1, double Q2, int order)
533{
534 //Define which stepping function should be used
535 const gsl_odeiv2_step_type * T = gsl_odeiv2_step_rk4;
536
537 //Allocate space for the stepping function
538 gsl_odeiv2_step * s = gsl_odeiv2_step_alloc(T, NumberOfRGEs);
539
540 //Define the absolute (A) and relative (R) error on y at each step.
541 //The real error will be compared to the following error estimate:
542 // A + R * |y_i|
543 gsl_odeiv2_control * c = gsl_odeiv2_control_y_new(1e-6, 0.0);
544
545 //Allocate space for the evolutor
546 gsl_odeiv2_evolve * e = gsl_odeiv2_evolve_alloc(NumberOfRGEs);
547
548 //Definition of the RGE system (the Jacobian is not necessary for the RK4 method; it's an empty function here)
549 gsl_odeiv2_system RGEsystem = {RGEsZ2, JacobianZ2, NumberOfRGEs, &order};
550
551 //Set starting and end point as natural logarithmic scales (conversion from decadic log scale)
552 double t1 = Q1*log(10.0);
553 double t2 = Q2*log(10.0);
554 double tNLOuni = NLOuniscale*log(10.0);
555
556 //Set initial step size
557 double InitialStepSize = 1e-6;
558
559 //Run!
560 while (t1 < t2)
561 {
562 int status = gsl_odeiv2_evolve_apply (e, c, s, &RGEsystem, &t1, t2, &InitialStepSize, InitialValues);
563 if(status != GSL_SUCCESS) break;
564
565 //intermediate checks if appropriate
566 if(RGEcheckZ2(InitialValues,t1,Rpeps,tNLOuni) != 0) break;
567 }
568
569 gsl_odeiv2_evolve_free (e);
570 gsl_odeiv2_control_free (c);
571 gsl_odeiv2_step_free (s);
572
573 //Return the decadic log scale at which the evolution stopped
574 return t1/log(10.0);
575}
576
578{
579 std::string modelflag = myGTHDMZ2->getZ2ModelType();
580 std::string RGEorder = myGTHDMZ2->getRGEorderflag();
581 int flag;
582 //flag will be used to transport information about model and RGEorder to the Runner:
583 //flag=3*(0 for type I, 1 for type II, 2 for type X, 3 for type Y, 4 for inert) + (0 for LO, 1 for approxNLO and 2 for NLO)
584 if( RGEorder == "LO" )
585 flag = 0;
586 else if( RGEorder == "approxNLO" )
587 flag = 1;
588 else if( RGEorder == "NLO" )
589 flag = 2;
590 else {
591 throw std::runtime_error("GeneralTHDMZ2Runner: RGEorder can be only any of \"LO\", \"approxNLO\" or \"NLO\"");
592 }
593
594 double g1_at_MZ = sqrt(4.0*M_PI*Ale/cW2);
595 double g2_at_MZ = sqrt(4.0*M_PI*Ale/(1-cW2));
596 double g3_at_MZ = sqrt(4.0*M_PI*Als);
597
598 double Ytop_at_MZ = (sqrt(2.0)*myGTHDMZ2->getQuarks(QCD::TOP).getMass())/(vev*sinb);
599 double Ybottom1_at_MZ = 0.0;
600 double Ybottom2_at_MZ = 0.0;
601 double Ytau1_at_MZ = 0.0;
602 double Ytau2_at_MZ = 0.0;
603
604 /*link these to the SM values*/
605 double Mb_at_MZ = 2.96;//GeV
606 double Mtau_at_MZ = 1.75;//GeV
607
608 if( modelflag == "type1" ) {
609 Ybottom2_at_MZ = (sqrt(2.0)*Mb_at_MZ)/(vev*sinb);
610 Ytau2_at_MZ = (sqrt(2.0)*Mtau_at_MZ)/(vev*sinb);
611 }
612 else if( modelflag == "type2" ) {
613 Ybottom1_at_MZ = (sqrt(2.0)*Mb_at_MZ)/(vev*cosb);
614 Ytau1_at_MZ = (sqrt(2.0)*Mtau_at_MZ)/(vev*cosb);
615
616 flag += 3;
617 }
618 else if( modelflag == "typeX" ) {
619 Ybottom2_at_MZ = (sqrt(2.0)*Mb_at_MZ)/(vev*sinb);
620 Ytau1_at_MZ = (sqrt(2.0)*Mtau_at_MZ)/(vev*cosb);
621
622 flag += 6;
623 }
624 else if( modelflag == "typeY" ) {
625 Ybottom1_at_MZ = (sqrt(2.0)*Mb_at_MZ)/(vev*cosb);
626 Ytau2_at_MZ = (sqrt(2.0)*Mtau_at_MZ)/(vev*sinb);
627 flag += 9;
628 }
629 else if( modelflag == "inert" ) {
630 Ytop_at_MZ = 0.0;
631
632 flag += 12;
633 }
634 else {
635 throw std::runtime_error("modelflag can be only any of \"type1\", \"type2\", \"typeX\", \"typeY\", \"inert\"");
636 }
637
638 double lambda1_at_MZ = myGTHDMZ2->getlambda1_Z2();
639 double lambda2_at_MZ = myGTHDMZ2->getlambda2_Z2();
640 double lambda3_at_MZ = myGTHDMZ2->getlambda3_Z2();
641 double lambda4_at_MZ = myGTHDMZ2->getlambda4_Z2();
642 double lambda5_at_MZ = myGTHDMZ2->getlambda5_Z2();
643
644 double lambda345 = lambda3_at_MZ + lambda4_at_MZ + lambda5_at_MZ;
645
646 double m12_2_at_MZ = m12_2;
647 double m11_2_at_MZ = tanb*m12_2_at_MZ - vev*vev*(cosb*cosb*lambda1_at_MZ - sinb*sinb*lambda345)/2.;
648 double m22_2_at_MZ = m12_2_at_MZ/tanb - vev*vev*(sinb*sinb*lambda2_at_MZ - cosb*cosb*lambda345)/2.;
649
650 if(fabs(Q_GTHDM-log10(MZ))<0.005) //at MZ scale
651 {
652 Q_cutoff = log10(MZ);
653
654 g1_at_Q = g1_at_MZ;
655 g2_at_Q = g2_at_MZ;
656 g3_at_Q = g3_at_MZ;
657 Ytop_at_Q = Ytop_at_MZ;
658 Ybottom1_at_Q = Ybottom1_at_MZ;
659 Ybottom2_at_Q = Ybottom2_at_MZ;
660 Ytau1_at_Q = Ytau1_at_MZ;
661 Ytau2_at_Q = Ytau2_at_MZ;
662 m11_2_at_Q = m11_2_at_MZ;
663 m22_2_at_Q = m22_2_at_MZ;
664 m12_2_at_Q = m12_2_at_MZ;
665 lambda1_at_Q = lambda1_at_MZ;
666 lambda2_at_Q = lambda2_at_MZ;
667 lambda3_at_Q = lambda3_at_MZ;
668 lambda4_at_Q = lambda4_at_MZ;
669 lambda5_at_Q = lambda5_at_MZ;
670 }
671 else //at some other scale
672 {
673 double InitVals[14];
674 InitVals[0] = g1_at_MZ;
675 InitVals[1] = g2_at_MZ;
676 InitVals[2] = g3_at_MZ;
677 InitVals[3] = Ytop_at_MZ;
678 InitVals[4] = Ybottom1_at_MZ+Ybottom2_at_MZ;
679 InitVals[5] = Ytau1_at_MZ+Ytau2_at_MZ;
680 InitVals[6] = m11_2_at_MZ;
681 InitVals[7] = m22_2_at_MZ;
682 InitVals[8] = m12_2_at_MZ;
683 InitVals[9] = lambda1_at_MZ;
684 InitVals[10] = lambda2_at_MZ;
685 InitVals[11] = lambda3_at_MZ;
686 InitVals[12] = lambda4_at_MZ;
687 InitVals[13] = lambda5_at_MZ;
688
689 // Running up to Q_cutoff <= Q_GTHDM
690 Q_cutoff = RGEGeneralTHDMZ2Runner(InitVals, 14, log10(MZ), Q_GTHDM, flag);
691
692 g1_at_Q = InitVals[0];
693 g2_at_Q = InitVals[1];
694 g3_at_Q = InitVals[2];
695 Ytop_at_Q = InitVals[3];
696 Ybottom1_at_Q = 0.0;
697 Ybottom2_at_Q = 0.0;
698 Ytau1_at_Q = 0.0;
699 Ytau2_at_Q = 0.0;
700
701 if( modelflag == "type1" ) {
702 Ybottom2_at_Q = InitVals[4];
703 Ytau2_at_Q = InitVals[5];
704 }
705 else if( modelflag == "type2" ) {
706 Ybottom1_at_Q = InitVals[4];
707 Ytau1_at_Q = InitVals[5];
708 }
709 else if( modelflag == "typeX" ) {
710 Ybottom2_at_Q = InitVals[4];
711 Ytau1_at_Q = InitVals[5];
712 }
713 else if( modelflag == "typeY" ) {
714 Ybottom1_at_Q = InitVals[4];
715 Ytau2_at_Q = InitVals[5];
716 }
717 else if( modelflag == "inert" ) {
718 Ytop_at_Q = 0.0;
719 }
720 else {
721 throw std::runtime_error("modelflag can be only any of \"type1\", \"type2\", \"typeX\" or \"typeY\"");
722 }
723
724 m11_2_at_Q = InitVals[6];
725 m22_2_at_Q = InitVals[7];
726 m12_2_at_Q = InitVals[8];
727 lambda1_at_Q = InitVals[9];
728 lambda2_at_Q = InitVals[10];
729 lambda3_at_Q = InitVals[11];
730 lambda4_at_Q = InitVals[12];
731 lambda5_at_Q = InitVals[13];
732 }
733}
734
735
736// Wavefunction renormalization contributions to unitarity, from Grinstein:2015rtl
738{
739 double WFRcomb1a;
740 double WFRcomb1b;
741 double WFRcomb2a;
742 double WFRcomb3a;
743 double WFRcomb3b;
744 double WFRcomb4a;
745
746 double B000mh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_0_0_mHh2(MZ2,mHl2).real();
747 double B000mH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_0_0_mHh2(MZ2,mHh2).real();
748 double B00mHpmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_0_mHp2_mHl2(MZ2,mHp2,mHl2).real();
749 double B00mHpmH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_0_mHp2_mHh2(MZ2,mHp2,mHh2).real();
750 double B00mAmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_0_mA2_mHl2(MZ2,mA2,mHl2).real();
751 double B00mAmH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_0_mA2_mHh2(MZ2,mA2,mHh2).real();
752 double B0mh00 = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_0_0(MZ2,mHl2).real();
753 double B0mh0mHp = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_0_mHp2(MZ2,mHl2,mHp2).real();
754 double B0mh0mA = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_0_mA2(MZ2,mHl2,mA2).real();
755 double B0mhmhmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_mHl2_mHl2(MZ2,mHl2).real();
756 double B0mhmHmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_mHh2_mHl2(MZ2,mHl2,mHh2).real();
757 double B0mhmHmH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_mHh2_mHh2(MZ2,mHl2,mHh2).real();
758 double B0mhmHpmHp = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_mHp2_mHp2(MZ2,mHl2,mHp2).real();
759 double B0mhmAmA = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHl2_mA2_mA2(MZ2,mHl2,mA2).real();
760 double B0mH00 = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_0_0(MZ2,mHh2).real();
761 double B0mH0mHp = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_0_mHp2(MZ2,mHh2,mHp2).real();
762 double B0mH0mA = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_0_mA2(MZ2,mHh2,mA2).real();
763 double B0mHmhmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_mHl2_mHl2(MZ2,mHh2,mHl2).real();
764 double B0mHmHmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_mHh2_mHl2(MZ2,mHh2,mHl2).real();
765 double B0mHmHmH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_mHh2_mHh2(MZ2,mHh2).real();
766 double B0mHmHpmHp = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_mHp2_mHp2(MZ2,mHh2,mHp2).real();
767 double B0mHmAmA = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHh2_mA2_mA2(MZ2,mHh2,mA2).real();
768 double B0mHp0mh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHp2_0_mHl2(MZ2,mHp2,mHl2).real();
769 double B0mHp0mH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHp2_0_mHh2(MZ2,mHp2,mHh2).real();
770 double B0mHpmHpmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHp2_mHp2_mHl2(MZ2,mHp2,mHl2).real();
771 double B0mHpmHpmH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mHp2_mHp2_mHh2(MZ2,mHp2,mHh2).real();
772 double B0mA0mh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mA2_0_mHl2(MZ2,mA2,mHl2).real();
773 double B0mA0mH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mA2_0_mHh2(MZ2,mA2,mHh2).real();
774 double B0mAmAmh = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mA2_mA2_mHl2(MZ2,mA2,mHl2).real();
775 double B0mAmAmH = myGTHDMZ2->getMyGTHDMCache()->B0_MZ2_mA2_mA2_mHh2(MZ2,mA2,mHh2).real();
776
777 double ddpB000mh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_0_mHl2(MZ2,mHl2).real();
778 double ddpB000mH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_0_mHh2(MZ2,mHh2).real();
779 double ddpB00mHpmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_mHp2_mHl2(MZ2,mHp2,mHl2).real();
780 double ddpB00mHpmH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_mHp2_mHh2(MZ2,mHp2,mHh2).real();
781 double ddpB00mHpmA = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_mHp2_mA2(MZ2,mHp2,mA2).real();
782 double ddpB00mAmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_mA2_mHl2(MZ2,mA2,mHl2).real();
783 double ddpB00mAmH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_0_mA2_mHh2(MZ2,mA2,mHh2).real();
784 double ddpB0mh00 = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_0_0(MZ2,mHl2).real();
785 double ddpB0mh0mHp = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_0_mHp2(MZ2,mHl2,mHp2).real();
786 double ddpB0mh0mA = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_0_mA2(MZ2,mHl2,mA2).real();
787 double ddpB0mhmhmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_mHl2_mHl2(MZ2,mHl2).real();
788 double ddpB0mhmHmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_mHh2_mHl2(MZ2,mHl2,mHh2).real();
789 double ddpB0mhmHmH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_mHh2_mHh2(MZ2,mHl2,mHh2).real();
790 double ddpB0mhmHpmHp = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_mHp2_mHp2(MZ2,mHl2,mHp2).real();
791 double ddpB0mhmAmA = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHl2_mA2_mA2(MZ2,mHl2,mA2).real();
792 double ddpB0mH00 = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_0_0(MZ2,mHh2).real();
793 double ddpB0mH0mHp = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_0_mHp2(MZ2,mHh2,mHp2).real();
794 double ddpB0mH0mA = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_0_mA2(MZ2,mHh2,mA2).real();
795 double ddpB0mHmhmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_mHl2_mHl2(MZ2,mHh2,mHl2).real();
796 double ddpB0mHmHmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_mHh2_mHl2(MZ2,mHh2,mHl2).real();
797 double ddpB0mHmHmH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_mHh2_mHh2(MZ2,mHh2).real();
798 double ddpB0mHmHpmHp = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_mHp2_mHp2(MZ2,mHh2,mHp2).real();
799 double ddpB0mHmAmA = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHh2_mA2_mA2(MZ2,mHh2,mA2).real();
800 double ddpB0mHp0mh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHp2_0_mHl2(MZ2,mHp2,mHl2).real();
801 double ddpB0mHp0mH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHp2_0_mHh2(MZ2,mHp2,mHh2).real();
802 double ddpB0mHp0mA = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHp2_0_mA2(MZ2,mHp2,mA2).real();
803 double ddpB0mHpmHpmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHp2_mHp2_mHl2(MZ2,mHp2,mHl2).real();
804 double ddpB0mHpmHpmH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mHp2_mHp2_mHh2(MZ2,mHp2,mHh2).real();
805 double ddpB0mA0mh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mA2_0_mHl2(MZ2,mA2,mHl2).real();
806 double ddpB0mA0mH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mA2_0_mHh2(MZ2,mA2,mHh2).real();
807 double ddpB0mA0mHp = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mA2_0_mHp2(MZ2,mA2,mHp2).real();
808 double ddpB0mAmAmh = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mA2_mA2_mHl2(MZ2,mA2,mHl2).real();
809 double ddpB0mAmAmH = myGTHDMZ2->getMyGTHDMCache()->B0p_MZ2_mA2_mA2_mHh2(MZ2,mA2,mHh2).real();
810
811 WFRcomb1a = 3.0*mHl2*mHl2*cosb*cosb*sin(bma)*sin(bma) * ddpB000mh
812 + 3.0*mHh2*mHh2*cos(bma)*cos(bma)*cosb*cosb * ddpB000mH
813 + 2.0*(mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma)*cosb*cosb * ddpB00mHpmh
814 + 2.0*(mHh2-mHp2)*(mHh2-mHp2)*sin(bma)*sin(bma)*cosb*cosb * ddpB00mHpmH
815 + 2.0*(mA2-mHp2)*(mA2-mHp2)*cosb*cosb * ddpB00mHpmA
816 + (mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma)*cosb*cosb * ddpB00mAmh
817 + (mA2-mHh2)*(mA2-mHh2)*cosb*cosb*sin(bma)*sin(bma) * ddpB00mAmH
818 + 1.5*mHl2*mHl2*sin(alpha)*sin(alpha)*sin(bma)*sin(bma) * ddpB0mh00
819 + 2.0*(mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma)*sin(alpha)*sin(alpha) * ddpB0mh0mHp
820 + (mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma)*sin(alpha)*sin(alpha) * ddpB0mh0mA
821 + 9.0*sin(alpha)*sin(alpha)*pow(-mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta))
822 +16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta),2)/(512.0*pow(cosb*sinb,4)) * ddpB0mhmhmh
823 + sin(alpha)*sin(alpha)*pow((cos(alpha)/sinb+sin(alpha)/cosb)*(m12_2+cos(alpha)*sin(alpha)*(mHh2+2.0*mHl2-(3.0*m12_2)/(cosb*sinb))),2) * ddpB0mhmHmh
824 + sin(alpha)*sin(alpha)*sin(bma)*sin(bma)*pow(-2.0*m12_2+(2.0*mHh2+mHl2-(3.0*m12_2)/(cosb*sinb))*sin(2.0*alpha),2)/(8.0*cosb*cosb*sinb*sinb) * ddpB0mhmHmH
825 + (sin(alpha)*sin(alpha)*pow((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
826 +cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb),2))/(64.0*pow(cosb*sinb,4)) * ddpB0mhmHpmHp
827 + (sin(alpha)*sin(alpha)*pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(8.0*m12_2-(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb),2))/(128.0*pow(cosb*sinb,4)) * ddpB0mhmAmA
828 + 1.5*mHh2*mHh2*cos(alpha)*cos(alpha)*cos(bma)*cos(bma) * ddpB0mH00
829 + 2.0*(mHh2-mHp2)*(mHh2-mHp2)*cos(alpha)*cos(alpha)*sin(bma)*sin(bma) * ddpB0mH0mHp
830 + (mA2-mHh2)*(mA2-mHh2)*cos(alpha)*cos(alpha)*sin(bma)*sin(bma) * ddpB0mH0mA
831 + cos(alpha)*cos(alpha)*cos(bma)*cos(bma)*pow(m12_2+cos(alpha)*sin(alpha)*(mHh2-3.0*m12_2/(cosb*sinb))+mHl2*sin(2.0*alpha),2)/(2.0*cosb*cosb*sinb*sinb) * ddpB0mHmhmh
832 + cos(alpha)*cos(alpha)*sin(bma)*sin(bma)*pow(m12_2*cosb*sinb+0.5*sin(2.0*alpha)*(3.0*m12_2-(2.0*mHh2+mHl2)*cosb*sinb),2)/pow(cosb*sinb,4) * ddpB0mHmHmh
833 + 9.0*cos(alpha)*cos(alpha)*pow(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
834 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta),2)/(512.0*pow(cosb*sinb,4)) * ddpB0mHmHmH
835 + (cos(alpha)*cos(alpha)*pow((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)
836 +2.0*(mHh2+2.0*mHp2)*cos(bma)-(3.0*mHh2+2.0*mHp2)*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta),2))/(256.0*pow(cosb*sinb,4)) * ddpB0mHmHpmHp
837 + (cos(alpha)*cos(alpha)*pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)
838 -2.0*(2.0*mA2+mHh2)*cos(bma)+(2.0*mA2+3.0*mHh2)*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta),2))/(512.0*pow(cosb*sinb,4)) * ddpB0mHmAmA
839 + 2.0*(mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma)*sinb*sinb * ddpB0mHp0mh
840 + 2.0*(mHh2-mHp2)*(mHh2-mHp2)*sin(bma)*sin(bma)*sinb*sinb * ddpB0mHp0mH
841 + 2.0*(mA2-mHp2)*(mA2-mHp2)*sinb*sinb * ddpB0mHp0mA
842 + 2.0*pow((m12_2*cos(alpha+beta))/(sinb*sinb*cosb*cosb)-(mHl2*cos(bma)*cos(2.0*beta))/(cosb*sinb)-(mHl2+2.0*mHp2)*sin(bma),2)*sinb*sinb * ddpB0mHpmHpmh
843 + 2.0*pow(sinb*((mHh2+2.0*mHp2)*cos(bma)-mHh2*cos(2.0*beta)*sin(bma)/(cosb*sinb)-m12_2*sin(alpha+beta)/(sinb*sinb*cosb*cosb)),2) * ddpB0mHpmHpmH
844 + (mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma)*sinb*sinb * ddpB0mA0mh
845 + (mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma)*sinb*sinb * ddpB0mA0mH
846 + 2.0*(mA2-mHp2)*(mA2-mHp2)*sinb*sinb * ddpB0mA0mHp
847 + pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(8.0*m12_2-(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb),2)/(64.0*pow(cosb,4)*sinb*sinb) * ddpB0mAmAmh
848 + pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+2.0*mA2*cos(alpha+3.0*beta)
849 + 3.0*mHh2*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta),2)/(256.0*pow(cosb,4)*sinb*sinb) * ddpB0mAmAmH;
850
851 WFRcomb1b = (mHl2*(mA2*(2.0*mHl2-3.0*mHp2)+mHl2*mHp2)*sin(2.0*bma)*2.0*sinb*cosb)/(2.0*mA2*mHp2) * B000mh
852 - (mHh2*(mA2*(2.0*mHh2-3.0*mHp2)+mHh2*mHp2)*sin(2.0*bma)*2.0*sinb*cosb)/(2.0*mA2*mHp2) * B000mH
853 + ((mHl2-mHp2)*cos(bma)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
854 +cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb)))/(2.0*mHp2*sinb*cosb) * B00mHpmh
855 + ((mHp2-mHh2)*sin(bma)*((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
856 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta)))/(4.0*mHp2*sinb*cosb) * B00mHpmH
857 + ((mHl2-mA2)*cos(bma)*((mHl2-2.0*mA2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
858 +cos(alpha+beta)*(-8.0*m12_2+(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb)))/(4.0*mA2*sinb*cosb) * B00mAmh
859 + ((mHh2-mA2)*sin(bma)*((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+2.0*mA2*cos(alpha+3.0*beta)
860 +3.0*mHh2*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta)))/(8.0*mA2*sinb*cosb) * B00mAmH
861 + (3.0*mHh2*mHl2*sin(2.0*alpha)*sin(2.0*bma))/(4.0*(mHh2-mHl2)) * B0mh00
862 + ((mHl2-mHp2)*(mHp2-mHh2)*sin(2.0*alpha)*sin(2.0*bma))/(mHh2-mHl2) * B0mh0mHp
863 - ((mA2-mHh2)*(mA2-mHl2)*sin(2.0*alpha)*sin(2.0*bma))/(2.0*(mHh2-mHl2)) * B0mh0mA
864 + 3.0*cos(bma)*sin(2.0*alpha)*(m12_2+cos(alpha)*sin(alpha)*(mHh2-(3.0*m12_2)/(sinb*cosb))+mHl2*sin(2.0*alpha))
865 *(-mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta))
866 +16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta))/(32.0*(mHl2-mHh2)*pow(sinb*cosb,3)) * B0mhmhmh
867 + (sin(2.0*bma)*sin(2.0*alpha)
868 *(4.0*m12_2*m12_2+2.0*m12_2*(mHl2-mHh2)*sin(2.0*alpha)
869 -(2.0*mHh2*mHh2+5.0*mHh2*mHl2+2.0*mHl2*mHl2-(9.0*m12_2*(mHh2+mHl2))/(sinb*cosb)
870 +(9.0*m12_2*m12_2)/(sinb*sinb*cosb*cosb))*sin(2.0*alpha)*sin(2.0*alpha)))/(8.0*(mHh2-mHl2)*sinb*sinb*cosb*cosb) * B0mhmHmh
871 - 3.0*sin(2.0*alpha)*sin(bma)*(m12_2-cos(alpha)*sin(alpha)*(2.0*mHh2+mHl2-(3.0*m12_2)/(sinb*cosb)))
872 *(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
873 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta))/(32.0*(mHh2-mHl2)*pow(sinb*cosb,3)) * B0mhmHmH
874 + sin(2.0*alpha)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb))
875 *((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
876 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(128.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mhmHpmHp
877 + sin(2.0*alpha)*((mHl2-2.0*mA2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb))
878 *((mHh2-2.0*mA2)*cos(alpha-5.0*beta)+2.0*(2.0*mA2+mHh2)*cos(bma)-2.0*mA2*cos(alpha+3.0*beta)
879 -3.0*mHh2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(256.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mhmAmA
880 - (3.0*mHh2*mHl2*sin(2.0*alpha)*sin(2.0*bma))/(4.0*(mHh2-mHl2)) * B0mH00
881 + ((mHh2-mHp2)*(mHl2-mHp2)*sin(2.0*alpha)*sin(2.0*bma))/(mHh2-mHl2) * B0mH0mHp
882 + ((mA2-mHh2)*(mA2-mHl2)*sin(2.0*alpha)*sin(2.0*bma))/(2.0*(mHh2-mHl2)) * B0mH0mA
883 + 3.0*cos(bma)*sin(2.0*alpha)*(m12_2+cos(alpha)*sin(alpha)*(mHh2-(3.0*m12_2)/(sinb*cosb))+mHl2*sin(2.0*alpha))
884 *(-mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta))
885 +16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta))/(32.0*(mHh2-mHl2)*pow(sinb*cosb,3)) * B0mHmhmh
886 + (sin(2.0*bma)*sin(2.0*alpha)*(-4.0*m12_2*m12_2+2.0*m12_2*(mHh2-mHl2)*sin(2.0*alpha)
887 +(2.0*mHh2*mHh2+5.0*mHh2*mHl2+2.0*mHl2*mHl2
888 -(9.0*m12_2*(mHh2+mHl2))/(sinb*cosb)
889 +(9.0*m12_2*m12_2)/(sinb*sinb*cosb*cosb))*sin(2.0*alpha)*sin(2.0*alpha)))
890 /(8.0*(mHh2-mHl2)*sinb*sinb*cosb*cosb) * B0mHmHmh
891 + 3.0*sin(bma)*sin(2.0*alpha)*(m12_2-cos(alpha)*sin(alpha)*(2.0*mHh2+mHl2-(3.0*m12_2)/(sinb*cosb)))
892 *(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
893 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta))/(32.0*(mHh2-mHl2)*pow(sinb*cosb,3)) * B0mHmHmH
894 - sin(2.0*alpha)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb))
895 *((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
896 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(128.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mHmHpmHp
897 - sin(2.0*alpha)*((mHl2-2.0*mA2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mA2)*2.0*sinb*cosb))
898 *((mHh2-2.0*mA2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mA2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
899 -2.0*mA2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(256.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mHmAmA
900 - (mHl2*(mHl2-mHp2)*sin(2.0*bma)*2.0*sinb*cosb)/mHp2 * B0mHp0mh
901 + (mHh2*(mHh2-mHp2)*sin(2.0*bma)*2.0*sinb*cosb)/mHp2 * B0mHp0mH
902 + ((mHp2-mHl2)*cos(bma)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
903 +cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb)))/(2.0*mHp2*sinb*cosb) * B0mHpmHpmh
904 + ((mHh2-mHp2)*sin(bma)*((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
905 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta)))/(4.0*mHp2*sinb*cosb) * B0mHpmHpmH
906 + (mHl2*(mA2-mHl2)*sin(bma)*cos(bma)*2.0*sinb*cosb)/mA2 * B0mA0mh
907 + ((mHh2-mA2)*mHh2*sin(bma)*cos(bma)*2.0*sinb*cosb)/mA2 * B0mA0mH
908 + ((mA2-mHl2)*cos(bma)*((-2.0*mA2+mHl2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
909 +cos(alpha+beta)*(-8.0*m12_2+(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb)))/(4.0*mA2*sinb*cosb) * B0mAmAmh
910 +((mA2-mHh2)*sin(bma)*((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+2.0*mA2*cos(alpha+3.0*beta)
911 +3.0*mHh2*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta)))/(8.0*mA2*sinb*cosb) * B0mAmAmH;
912
913 WFRcomb2a = 1.5*mHl2*mHl2*sin(bma)*sin(bma) * ddpB000mh
914 + 1.5*mHh2*mHh2*cos(bma)*cos(bma) * ddpB000mH
915 + (mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma) * ddpB00mHpmh
916 + (mHh2-mHp2)*(mHh2-mHp2)*sin(bma)*sin(bma) * ddpB00mHpmH
917 + (mA2-mHp2)*(mA2-mHp2) * ddpB00mHpmA
918 + 0.5*(mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma) * ddpB00mAmh
919 + 0.5*(mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma) * ddpB00mAmH
920 + 0.75*mHl2*mHl2*sin(bma)*sin(bma) * ddpB0mh00
921 + (mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma) * ddpB0mh0mHp
922 + 0.5*(mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma) * ddpB0mh0mA
923 + 9.0*pow(16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta)
924 -mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta)),2)/(1024.0*pow(cosb*sinb,4)) * ddpB0mhmhmh
925 + 0.5*pow(cos(alpha)/sinb + sin(alpha)/cosb,2)
926 *pow(m12_2+cos(alpha)*sin(alpha)*(mHh2+2.0*mHl2-3.0*m12_2/(cosb*sinb)),2) * ddpB0mhmHmh
927 + (pow(-2.0*m12_2+(2.0*mHh2+mHl2-3.0*m12_2/(cosb*sinb))*sin(2.0*alpha),2)
928 *sin(bma)*sin(bma))/(16.0*cosb*cosb*sinb*sinb) * ddpB0mhmHmH
929 + pow((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*cosb*sinb
930 +cos(alpha+beta)*(-4.0*m12_2+(3.0*mHl2+2.0*mHp2)*cosb*sinb),2)/(32.0*pow(cosb*sinb,4)) * ddpB0mhmHpmHp
931 + pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*cosb*sinb
932 +cos(alpha+beta)*(4.0*m12_2-(2.0*mA2+3.0*mHl2)*cosb*sinb),2)/(64.0*pow(cosb*sinb,4)) * ddpB0mhmAmA
933 + 0.75*mHh2*mHh2*cos(bma)*cos(bma) * ddpB0mH00
934 + (mHh2 - mHp2)*(mHh2 - mHp2)*sin(bma)*sin(bma) * ddpB0mH0mHp
935 + 0.5*(mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma) * ddpB0mH0mA
936 + 0.25*pow(cos(alpha)/sinb + sin(alpha)/cosb,2)
937 *pow(m12_2+cos(alpha)*sin(alpha)*(mHh2+2.0*mHl2-3.0*m12_2/(cosb*sinb)),2) * ddpB0mHmhmh
938 + (pow(-2.0*m12_2+(2.0*mHh2+mHl2-3.0*m12_2/(cosb*sinb))*sin(2.0*alpha),2)
939 *sin(bma)*sin(bma))/(8.0*cosb*cosb*sinb*sinb) * ddpB0mHmHmh
940 + 9.0*pow(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
941 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta),2)/(1024.0*pow(cosb*sinb,4)) * ddpB0mHmHmH
942 + pow((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-(3.0*mHh2+2.0*mHp2)*cos(alpha+3.0*beta)
943 -16.0*m12_2*sin(alpha+beta),2)/(512*pow(cosb*sinb,4)) * ddpB0mHmHpmHp
944 + pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+(2.0*mA2+3.0*mHh2)*cos(alpha+3.0*beta)
945 +16.0*m12_2*sin(alpha+beta),2)/(1024*pow(cosb*sinb,4)) * ddpB0mHmAmA
946 + (mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma) * ddpB0mHp0mh
947 + (mHh2-mHp2)*(mHh2-mHp2)*sin(bma)*sin(bma) * ddpB0mHp0mH
948 + (mA2-mHp2)*(mA2-mHp2) * ddpB0mHp0mA
949 + pow((m12_2*cos(alpha+beta))/(cosb*cosb*sinb*sinb)-(mHl2*cos(bma)*cos(2.0*beta))/(cosb*sinb)
950 -(mHl2+2.0*mHp2)*sin(bma),2) * ddpB0mHpmHpmh
951 + pow((mHh2+2.0*mHp2)*cos(bma)-(mHh2*cos(2.0*beta)*sin(bma))/(cosb*sinb)
952 -(m12_2*sin(alpha+beta))/(cosb*cosb*sinb*sinb),2) * ddpB0mHpmHpmH
953 + 0.5*(mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma) * ddpB0mA0mh
954 + 0.5*(mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma) * ddpB0mA0mH
955 + (mA2-mHp2)*(mA2-mHp2) * ddpB0mA0mHp
956 + pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*2.0*cosb*sinb
957 +cos(alpha+beta)*(8.0*m12_2-(2.0*mA2+3.0*mHl2)*2.0*cosb*sinb),2)/(128.0*pow(cosb*sinb,4)) * ddpB0mAmAmh
958 + pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+(2.0*mA2+3.0*mHh2)*cos(alpha+3.0*beta)
959 +16.0*m12_2*sin(alpha+beta),2)/(512.0*pow(cosb*sinb,4)) * ddpB0mAmAmH;
960
961 WFRcomb3a = 0.5*mHl2*mHl2*(3.0-cos(2.0*beta))*sin(bma)*sin(bma) * ddpB000mh
962 + 0.5*mHh2*mHh2*(3.0-cos(2.0*beta))*cos(bma)*cos(bma) * ddpB000mH
963 + 2.0*(mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma)*sinb*sinb * ddpB00mHpmh
964 + 2.0*(mHh2-mHp2)*(mHh2-mHp2)*sin(bma)*sin(bma)*sinb*sinb * ddpB00mHpmH
965 + 2.0*(mA2-mHp2)*(mA2-mHp2)*sinb*sinb * ddpB00mHpmA
966 + (mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma)*cosb*cosb * ddpB00mAmh
967 + (mA2-mHh2)*(mA2-mHh2)*cosb*cosb*sin(bma)*sin(bma) * ddpB00mAmH
968 + 1.5*mHl2*mHl2*sin(alpha)*sin(alpha)*sin(bma)*sin(bma) * ddpB0mh00
969 + 2.0*(mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma)*sin(alpha)*sin(alpha) * ddpB0mh0mHp
970 + (mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma)*sin(alpha)*sin(alpha) * ddpB0mh0mA
971 + 9.0*sin(alpha)*sin(alpha)*pow(-mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta))
972 +16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta),2)/(512.0*pow(cosb*sinb,4)) * ddpB0mhmhmh
973 + sin(alpha)*sin(alpha)*pow((cos(alpha)/sinb+sin(alpha)/cosb)*(m12_2+cos(alpha)*sin(alpha)*(mHh2+2.0*mHl2-(3.0*m12_2)/(cosb*sinb))),2) * ddpB0mhmHmh
974 + sin(alpha)*sin(alpha)*sin(bma)*sin(bma)*pow(-2.0*m12_2+(2.0*mHh2+mHl2-(3.0*m12_2)/(cosb*sinb))*sin(2.0*alpha),2)/(8.0*cosb*cosb*sinb*sinb) * ddpB0mhmHmH
975 + (sin(alpha)*sin(alpha)*pow((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
976 +cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb),2))/(64.0*pow(cosb*sinb,4)) * ddpB0mhmHpmHp
977 + (sin(alpha)*sin(alpha)*pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(8.0*m12_2-(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb),2))/(128.0*pow(cosb*sinb,4)) * ddpB0mhmAmA
978 + 1.5*mHh2*mHh2*cos(alpha)*cos(alpha)*cos(bma)*cos(bma) * ddpB0mH00
979 + 2.0*(mHh2-mHp2)*(mHh2-mHp2)*cos(alpha)*cos(alpha)*sin(bma)*sin(bma) * ddpB0mH0mHp
980 + (mA2-mHh2)*(mA2-mHh2)*cos(alpha)*cos(alpha)*sin(bma)*sin(bma) * ddpB0mH0mA
981 + cos(alpha)*cos(alpha)*cos(bma)*cos(bma)*pow(m12_2+cos(alpha)*sin(alpha)*(mHh2-3.0*m12_2/(cosb*sinb))+mHl2*sin(2.0*alpha),2)/(2.0*cosb*cosb*sinb*sinb) * ddpB0mHmhmh
982 + cos(alpha)*cos(alpha)*sin(bma)*sin(bma)*pow(m12_2*cosb*sinb+0.5*sin(2.0*alpha)*(3.0*m12_2-(2.0*mHh2+mHl2)*cosb*sinb),2)/pow(cosb*sinb,4) * ddpB0mHmHmh
983 + 9.0*cos(alpha)*cos(alpha)*pow(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
984 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta),2)/(512.0*pow(cosb*sinb,4)) * ddpB0mHmHmH
985 + (cos(alpha)*cos(alpha)*pow((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)
986 +2.0*(mHh2+2.0*mHp2)*cos(bma)-(3.0*mHh2+2.0*mHp2)*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta),2))/(256.0*pow(cosb*sinb,4)) * ddpB0mHmHpmHp
987 + (cos(alpha)*cos(alpha)*pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)
988 -2.0*(2.0*mA2+mHh2)*cos(bma)+(2.0*mA2+3.0*mHh2)*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta),2))/(512.0*pow(cosb*sinb,4)) * ddpB0mHmAmA
989 + 2.0*(mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma)*cosb*cosb * ddpB0mHp0mh
990 + 2.0*(mHh2-mHp2)*(mHh2-mHp2)*sin(bma)*sin(bma)*cosb*cosb * ddpB0mHp0mH
991 + 2.0*(mA2-mHp2)*(mA2-mHp2)*cosb*cosb * ddpB0mHp0mA
992 + 2.0*cosb*cosb*pow((m12_2*cos(alpha+beta))/(cosb*cosb*sinb*sinb)-(mHl2*cos(bma)*cos(2.0*beta))/(cosb*sinb)-(mHl2+2.0*mHp2)*sin(bma),2) * ddpB0mHpmHpmh
993 + 2.0*cosb*cosb*pow((mHh2+2.0*mHp2)*cos(bma)-(mHh2*cos(2.0*beta)*sin(bma))/(cosb*sinb)-(m12_2*sin(alpha+beta))/(cosb*cosb*sinb*sinb),2) * ddpB0mHpmHpmH
994 + (mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma)*sinb*sinb * ddpB0mA0mh
995 + (mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma)*sinb*sinb * ddpB0mA0mH
996 + 2.0*(mA2-mHp2)*(mA2-mHp2)*sinb*sinb * ddpB0mA0mHp
997 + pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(8.0*m12_2-(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb),2)/(64.0*pow(cosb,4)*sinb*sinb) * ddpB0mAmAmh
998 + pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+2.0*mA2*cos(alpha+3.0*beta)
999 + 3.0*mHh2*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta),2)/(256.0*pow(cosb,4)*sinb*sinb) * ddpB0mAmAmH;
1000
1001 WFRcomb3b = ((mHl2*mHl2*mHp2+mA2*(-2.0*mHl2*mHl2+mHl2*mHp2))*sin(2.0*bma)*sinb*cosb)/(mA2*mHp2) * B000mh
1002 + ((-mHh2*mHh2*mHp2+mA2*(2.0*mHh2*mHh2-mHh2*mHp2))*sin(2.0*bma)*sinb*cosb)/(mA2*mHp2) * B000mH
1003 + ((mHp2-mHl2)*cos(bma)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
1004 +cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb)))/(2.0*mHp2*sinb*cosb) * B00mHpmh
1005 + ((mHh2-mHp2)*sin(bma)*((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
1006 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta)))/(4.0*mHp2*sinb*cosb) * B00mHpmH
1007 + ((mHl2-mA2)*cos(bma)*((mHl2-2.0*mA2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
1008 +cos(alpha+beta)*(-8.0*m12_2+(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb)))/(4.0*mA2*sinb*cosb) * B00mAmh
1009 + ((mHh2-mA2)*sin(bma)*((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+2.0*mA2*cos(alpha+3.0*beta)
1010 +3.0*mHh2*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta)))/(8.0*mA2*sinb*cosb) * B00mAmH
1011 + (3.0*mHh2*mHl2*sin(2.0*alpha)*sin(2.0*bma))/(4.0*(mHh2-mHl2)) * B0mh00
1012 + ((mHl2-mHp2)*(mHp2-mHh2)*sin(2.0*alpha)*sin(2.0*bma))/(mHh2-mHl2) * B0mh0mHp
1013 - ((mA2-mHh2)*(mA2-mHl2)*sin(2.0*alpha)*sin(2.0*bma))/(2.0*(mHh2-mHl2)) * B0mh0mA
1014 + 3.0*cos(bma)*sin(2.0*alpha)*(m12_2+cos(alpha)*sin(alpha)*(mHh2-(3.0*m12_2)/(sinb*cosb))+mHl2*sin(2.0*alpha))
1015 *(-mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta))
1016 +16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta))/(32.0*(mHl2-mHh2)*pow(sinb*cosb,3)) * B0mhmhmh
1017 + (sin(2.0*bma)*sin(2.0*alpha)
1018 *(4.0*m12_2*m12_2+2.0*m12_2*(mHl2-mHh2)*sin(2.0*alpha)
1019 -(2.0*mHh2*mHh2+5.0*mHh2*mHl2+2.0*mHl2*mHl2-(9.0*m12_2*(mHh2+mHl2))/(sinb*cosb)
1020 +(9.0*m12_2*m12_2)/(sinb*sinb*cosb*cosb))*sin(2.0*alpha)*sin(2.0*alpha)))/(8.0*(mHh2-mHl2)*sinb*sinb*cosb*cosb) * B0mhmHmh
1021 - 3.0*sin(2.0*alpha)*sin(bma)*(m12_2-cos(alpha)*sin(alpha)*(2.0*mHh2+mHl2-(3.0*m12_2)/(sinb*cosb)))
1022 *(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
1023 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta))/(32.0*(mHh2-mHl2)*pow(sinb*cosb,3)) * B0mhmHmH
1024 + sin(2.0*alpha)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb))
1025 *((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
1026 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(128.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mhmHpmHp
1027 + sin(2.0*alpha)*((mHl2-2.0*mA2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb))
1028 *((mHh2-2.0*mA2)*cos(alpha-5.0*beta)+2.0*(2.0*mA2+mHh2)*cos(bma)-2.0*mA2*cos(alpha+3.0*beta)
1029 -3.0*mHh2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(256.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mhmAmA
1030 - (3.0*mHh2*mHl2*sin(2.0*alpha)*sin(2.0*bma))/(4.0*(mHh2-mHl2)) * B0mH00
1031 + ((mHh2-mHp2)*(mHl2-mHp2)*sin(2.0*alpha)*sin(2.0*bma))/(mHh2-mHl2) * B0mH0mHp
1032 + ((mA2-mHh2)*(mA2-mHl2)*sin(2.0*alpha)*sin(2.0*bma))/(2.0*(mHh2-mHl2)) * B0mH0mA
1033 + 3.0*cos(bma)*sin(2.0*alpha)*(m12_2+cos(alpha)*sin(alpha)*(mHh2-(3.0*m12_2)/(sinb*cosb))+mHl2*sin(2.0*alpha))
1034 *(-mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta))
1035 +16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta))/(32.0*(mHh2-mHl2)*pow(sinb*cosb,3)) * B0mHmhmh
1036 + (sin(2.0*bma)*sin(2.0*alpha)*(-4.0*m12_2*m12_2+2.0*m12_2*(mHh2-mHl2)*sin(2.0*alpha)
1037 +(2.0*mHh2*mHh2+5.0*mHh2*mHl2+2.0*mHl2*mHl2
1038 -(9.0*m12_2*(mHh2+mHl2))/(sinb*cosb)
1039 +(9.0*m12_2*m12_2)/(sinb*sinb*cosb*cosb))*sin(2.0*alpha)*sin(2.0*alpha)))
1040 /(8.0*(mHh2-mHl2)*sinb*sinb*cosb*cosb) * B0mHmHmh
1041 + 3.0*sin(bma)*sin(2.0*alpha)*(m12_2-cos(alpha)*sin(alpha)*(2.0*mHh2+mHl2-(3.0*m12_2)/(sinb*cosb)))
1042 *(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
1043 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta))/(32.0*(mHh2-mHl2)*pow(sinb*cosb,3)) * B0mHmHmH
1044 - sin(2.0*alpha)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb))
1045 *((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
1046 -2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(128.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mHmHpmHp
1047 - sin(2.0*alpha)*((mHl2-2.0*mA2)*cos(alpha-3.0*beta)*2.0*sinb*cosb+cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mA2)*2.0*sinb*cosb))
1048 *((mHh2-2.0*mA2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mA2)*cos(bma)-3.0*mHh2*cos(alpha+3.0*beta)
1049 -2.0*mA2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(256.0*(mHh2-mHl2)*pow(cosb*sinb,4)) * B0mHmAmA
1050 + (mHl2*(mHl2-mHp2)*sin(2.0*bma)*2.0*sinb*cosb)/mHp2 * B0mHp0mh
1051 - (mHh2*(mHh2-mHp2)*sin(2.0*bma)*2.0*sinb*cosb)/mHp2 * B0mHp0mH
1052 + (mHl2-mHp2)*cos(bma)*((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
1053 +cos(alpha+beta)*(-8.0*m12_2+(3.0*mHl2+2.0*mHp2)*2.0*sinb*cosb))/(2.0*mHp2*sinb*cosb) * B0mHpmHpmh
1054 - (mHh2-mHp2)*sin(bma)*((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)
1055 -3.0*mHh2*cos(alpha+3.0*beta)-2.0*mHp2*cos(alpha+3.0*beta)-16.0*m12_2*sin(alpha+beta))/(4.0*mHp2*sinb*cosb) * B0mHpmHpmH
1056 + (mHl2*(mA2-mHl2)*sin(bma)*cos(bma)*2.0*sinb*cosb)/mA2 * B0mA0mh
1057 + ((mHh2-mA2)*mHh2*sin(bma)*cos(bma)*2.0*sinb*cosb)/mA2 * B0mA0mH
1058 + ((mA2-mHl2)*cos(bma)*((-2.0*mA2+mHl2)*cos(alpha-3.0*beta)*2.0*sinb*cosb
1059 +cos(alpha+beta)*(-8.0*m12_2+(2.0*mA2+3.0*mHl2)*2.0*sinb*cosb)))/(4.0*mA2*sinb*cosb) * B0mAmAmh
1060 +((mA2-mHh2)*sin(bma)*((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+2.0*mA2*cos(alpha+3.0*beta)
1061 +3.0*mHh2*cos(alpha+3.0*beta)+16.0*m12_2*sin(alpha+beta)))/(8.0*mA2*sinb*cosb) * B0mAmAmH;
1062
1063 WFRcomb4a = 0.5*mHl2*mHl2*sin(bma)*sin(bma) * ddpB000mh
1064 + 0.5*mHh2*mHh2*cos(bma)*cos(bma) * ddpB000mH
1065 + 0.5*(mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma) * ddpB00mAmh
1066 + 0.5*(mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma) * ddpB00mAmH
1067 + 0.75*mHl2*mHl2*sin(bma)*sin(bma) * ddpB0mh00
1068 + (mHl2-mHp2)*(mHl2-mHp2)*cos(bma)*cos(bma) * ddpB0mh0mHp
1069 + 0.5*(mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma) * ddpB0mh0mA
1070 + 9.0*pow(16.0*m12_2*cos(bma)*cos(bma)*cos(alpha+beta)
1071 -mHl2*(3.0*sin(bma)+sin(3.0*bma)+sin(3.0*alpha+beta)+3.0*sin(alpha+3.0*beta)),2)/(1024.0*pow(cosb*sinb,4)) * ddpB0mhmhmh
1072 + 0.5*pow(cos(alpha)/sinb + sin(alpha)/cosb,2)
1073 *pow(m12_2+cos(alpha)*sin(alpha)*(mHh2+2.0*mHl2-3.0*m12_2/(cosb*sinb)),2) * ddpB0mhmHmh
1074 + (pow(-2.0*m12_2+(2.0*mHh2+mHl2-3.0*m12_2/(cosb*sinb))*sin(2.0*alpha),2)
1075 *sin(bma)*sin(bma))/(16.0*cosb*cosb*sinb*sinb) * ddpB0mhmHmH
1076 + pow((mHl2-2.0*mHp2)*cos(alpha-3.0*beta)*cosb*sinb
1077 +cos(alpha+beta)*(-4.0*m12_2+(3.0*mHl2+2.0*mHp2)*cosb*sinb),2)/(32.0*pow(cosb*sinb,4)) * ddpB0mhmHpmHp
1078 + pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*cosb*sinb
1079 +cos(alpha+beta)*(4.0*m12_2-(2.0*mA2+3.0*mHl2)*cosb*sinb),2)/(64.0*pow(cosb*sinb,4)) * ddpB0mhmAmA
1080 + 0.75*mHh2*mHh2*cos(bma)*cos(bma) * ddpB0mH00
1081 + (mHh2 - mHp2)*(mHh2 - mHp2)*sin(bma)*sin(bma) * ddpB0mH0mHp
1082 + 0.5*(mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma) * ddpB0mH0mA
1083 + 0.25*pow(cos(alpha)/sinb + sin(alpha)/cosb,2)
1084 *pow(m12_2+cos(alpha)*sin(alpha)*(mHh2+2.0*mHl2-3.0*m12_2/(cosb*sinb)),2) * ddpB0mHmhmh
1085 + (pow(-2.0*m12_2+(2.0*mHh2+mHl2-3.0*m12_2/(cosb*sinb))*sin(2.0*alpha),2)
1086 *sin(bma)*sin(bma))/(8.0*cosb*cosb*sinb*sinb) * ddpB0mHmHmh
1087 + 9.0*pow(mHh2*(-3.0*cos(bma)+cos(3.0*bma)-cos(3.0*alpha+beta)+3.0*cos(alpha+3.0*beta))
1088 +16.0*m12_2*sin(bma)*sin(bma)*sin(alpha+beta),2)/(1024.0*pow(cosb*sinb,4)) * ddpB0mHmHmH
1089 + pow((mHh2-2.0*mHp2)*cos(alpha-5.0*beta)+2.0*(mHh2+2.0*mHp2)*cos(bma)-(3.0*mHh2+2.0*mHp2)*cos(alpha+3.0*beta)
1090 -16.0*m12_2*sin(alpha+beta),2)/(512*pow(cosb*sinb,4)) * ddpB0mHmHpmHp
1091 + pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+(2.0*mA2+3.0*mHh2)*cos(alpha+3.0*beta)
1092 +16.0*m12_2*sin(alpha+beta),2)/(1024*pow(cosb*sinb,4)) * ddpB0mHmAmA
1093 + 0.5*(mA2-mHl2)*(mA2-mHl2)*cos(bma)*cos(bma) * ddpB0mA0mh
1094 + 0.5*(mA2-mHh2)*(mA2-mHh2)*sin(bma)*sin(bma) * ddpB0mA0mH
1095 + (mA2-mHp2)*(mA2-mHp2) * ddpB0mA0mHp
1096 + pow((2.0*mA2-mHl2)*cos(alpha-3.0*beta)*2.0*cosb*sinb
1097 +cos(alpha+beta)*(8.0*m12_2-(2.0*mA2+3.0*mHl2)*2.0*cosb*sinb),2)/(128.0*pow(cosb*sinb,4)) * ddpB0mAmAmh
1098 + pow((2.0*mA2-mHh2)*cos(alpha-5.0*beta)-2.0*(2.0*mA2+mHh2)*cos(bma)+(2.0*mA2+3.0*mHh2)*cos(alpha+3.0*beta)
1099 +16.0*m12_2*sin(alpha+beta),2)/(512.0*pow(cosb*sinb,4)) * ddpB0mAmAmH;
1100
1101 WFRcomb1 = -(WFRcomb1a+WFRcomb1b)/(vev*vev);
1102 WFRcomb2 = -WFRcomb2a/(vev*vev);
1103 WFRcomb3 = -(WFRcomb3a+WFRcomb3b)/(vev*vev);
1104 WFRcomb4 = -WFRcomb4a/(vev*vev);
1105}
1106
1107
1109{
1110 vev = myGTHDMZ2->v();
1111 cW2 = myGTHDMZ2->c02();
1112 Ale = myGTHDMZ2->getAle();
1113 Als = myGTHDMZ2->getAlsMz();
1114 MZ = myGTHDMZ2->getMz();
1115 MZ2 = MZ*MZ;
1116
1118 mHl2 = myGTHDMZ2->getmH1sq();
1119 mHh2 = myGTHDMZ2->getmH2sq();
1120 mA2 = myGTHDMZ2->getmH3sq();
1121 mHp2 = myGTHDMZ2->getmHp2();
1122
1127 bma = myGTHDMZ2->getbma_Z2();
1128 alpha = beta-bma;
1129
1130 Q_GTHDM = myGTHDMZ2->getQ_GTHDM();
1131 Rpeps = myGTHDMZ2->getRpepsGTHDM();
1132 NLOuniscale = myGTHDMZ2->getNLOuniscaleGTHDM();
1133
1135
1136 // Wavefunction renormalization contributions to unitarity, from Grinstein:2015rtl
1137 WFRcomb1 = 0.;
1138 WFRcomb2 = 0.;
1139 WFRcomb3 = 0.;
1140 WFRcomb4 = 0.;
1141
1143 computeWFR_Z2();
1144
1145 // 1st row: (lambda1, lambda2, lambda3, lambda4, lambda5)
1146 // 2nd row: (Ytop, Ybottom1, Ybottom2, Ytau1, Ytau2)
1147 // 3rd row: (WFRcomb1, WFRcomb2, WFRcomb3, WFRcomb4, 0.0)
1148 myZ2_at_Q.assign(0, 0, lambda1_at_Q);
1149 myZ2_at_Q.assign(0, 1, lambda2_at_Q);
1150 myZ2_at_Q.assign(0, 2, lambda3_at_Q);
1151 myZ2_at_Q.assign(0, 3, lambda4_at_Q);
1152 myZ2_at_Q.assign(0, 4, lambda5_at_Q);
1153 myZ2_at_Q.assign(1, 0, Ytop_at_Q);
1154 myZ2_at_Q.assign(1, 1, Ybottom1_at_Q);
1155 myZ2_at_Q.assign(1, 2, Ybottom2_at_Q);
1156 myZ2_at_Q.assign(1, 3, Ytau1_at_Q);
1157 myZ2_at_Q.assign(1, 4, Ytau2_at_Q);
1158 myZ2_at_Q.assign(2, 0, WFRcomb1);
1159 myZ2_at_Q.assign(2, 1, WFRcomb2);
1160 myZ2_at_Q.assign(2, 2, WFRcomb3);
1161 myZ2_at_Q.assign(2, 3, WFRcomb4);
1162
1163 return (myZ2_at_Q);
1164}
int JacobianZ2(double t, const double y[], double *dfdy, double dfdt[], void *order)
int RGEsZ2(double t, const double y[], double beta[], void *flags)
int RGEcheckZ2(const double InitialValues[], const double t1, const double Rpeps, const double tNLOuni)
double getlambda2_Z2() const
A getter for the parameter of the scalar potential .
double getlambda4_Z2() const
A getter for the parameter of the scalar potential .
double getsinb_Z2() const
A getter for the sine of .
double getcosb_Z2() const
A getter for the cosine of .
double gettanb_Z2() const
A getter for the tangent of .
double getlambda5_Z2() const
A getter for the parameter of the scalar potential .
bool getWFRflag_Z2() const
A getter for the chosen option to include WFR in NLO unitarity.
double getlambda3_Z2() const
A getter for the parameter of the scalar potential .
double getbeta_Z2() const
A getter for .
Definition: GeneralTHDMZ2.h:91
double getlambda1_Z2() const
A getter for the parameter of the scalar potential .
std::string getZ2ModelType() const
A getter for the chosen Z2-symmetric model type.
double getbma_Z2() const
A getter for the difference between and .
double getm12sq_Z2() const
A getter for the parameter of the scalar potential .
virtual ~GeneralTHDMZ2Runner()
Runner destructor.
GeneralTHDMZ2Runner(const StandardModel &SM_i)
GeneralTHDMZ2Runner constructor.
gslpp::matrix< double > myZ2_at_Q
gslpp::matrix< double > getGTHDMZ2_at_Q()
The public function which contains all relevant GTHDMZ2 parameter after running.
virtual double RGEGeneralTHDMZ2Runner(double InitialValues[], unsigned long int NumberOfRGEs, double Q1, double Q2, int order)
const GeneralTHDMZ2 * myGTHDMZ2
@ TOP
Definition: QCD.h:328
A model class for the Standard Model.
An observable class for the quartic Higgs potential coupling combination .
An observable class for the quadratic Higgs potential coupling .
An observable class for the quadratic Higgs potential coupling .
Test Observable.
Test Observable.