a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
GeneralTHDMZ2Runner.cpp File Reference

Go to the source code of this file.

Functions

int JacobianZ2 (double t, const double y[], double *dfdy, double dfdt[], void *order)
 
int RGEcheckZ2 (const double InitialValues[], const double t1, const double Rpeps, const double tNLOuni)
 
int RGEsZ2 (double t, const double y[], double beta[], void *flags)
 

Function Documentation

◆ JacobianZ2()

int JacobianZ2 ( double  t,
const double  y[],
double *  dfdy,
double  dfdt[],
void *  order 
)

Definition at line 397 of file GeneralTHDMZ2Runner.cpp.

398{
399 return 0;
400}

◆ RGEcheckZ2()

int RGEcheckZ2 ( const double  InitialValues[],
const double  t1,
const double  Rpeps,
const double  tNLOuni 
)

Definition at line 402 of file GeneralTHDMZ2Runner.cpp.

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}

◆ RGEsZ2()

int RGEsZ2 ( double  t,
const double  y[],
double  beta[],
void *  flags 
)

Definition at line 21 of file GeneralTHDMZ2Runner.cpp.

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}
An observable class for the quadratic Higgs potential coupling .
An observable class for the quadratic Higgs potential coupling .
Test Observable.