1:Evaluate: BeginPackage["GM2Calc`"]
3:Evaluate: loopOrder::usage =
4 "Loop order of the calculation. Possible values: 0, 1, 2";
6:Evaluate: tanBetaResummation::usage =
7 "Enable/disable tan(beta) resummation. Possible values: True, False";
9:Evaluate: forceOutput::usage =
10 "Enforce output, even in case an error has occurred.";
12:Evaluate: runningCouplings::usage =
13 "Use running couplings in the THDM.";
15:Evaluate: alphaMZ::usage =
16 "Electromagnetic coupling in the MS-bar scheme at the scale MZ."
18:Evaluate: alpha0::usage =
19 "Electromagnetic coupling in the Thomson limit."
21:Evaluate: alphaS::usage =
31 "Top quark pole mass."
33:Evaluate: mcmc::usage =
34 "Charm quark MS-bar mass mc at the scale mc."
36:Evaluate: mu2GeV::usage =
37 "Up quark MS-bar mass mu at the scale 2 GeV."
39:Evaluate: mbmb::usage =
40 "Bottom quark MS-bar mass mb at the scale mb."
42:Evaluate: ms2GeV::usage =
43 "Strange quark MS-bar mass ms at the scale 2 GeV."
45:Evaluate: md2GeV::usage =
46 "Down quark MS-bar mass ms at the scale 2 GeV."
48:Evaluate: mbMZ::usage =
49 "Bottom quark DR-bar mass at the scale MZ."
52 "Tau lepton pole mass."
60:Evaluate: Mv1::usage =
61 "Lightest neutrino pole mass."
63:Evaluate: Mv2::usage =
64 "2nd lightest neutrino pole mass."
66:Evaluate: Mv3::usage =
67 "Heaviest neutrino pole mass."
69:Evaluate: CKM::usage =
72:Evaluate: amu::usage =
73 "Calculated value of the anomalous magnetic moment of the muon, amu = (g-2)/2.";
75:Evaluate: amu1L::usage =
76 "Calculated value of the anomalous magnetic moment of the muon, amu = (g-2)/2 (1-loop contribution only).";
78:Evaluate: amu2L::usage =
79 "Calculated value of the anomalous magnetic moment of the muon, amu = (g-2)/2 (2-loop contribution only).";
81:Evaluate: amu2LF::usage =
82 "Calculated value of the anomalous magnetic moment of the muon, amu = (g-2)/2 (2-loop fermionic contribution only).";
84:Evaluate: amu2LB::usage =
85 "Calculated value of the anomalous magnetic moment of the muon, amu = (g-2)/2 (2-loop bosonic contribution only).";
87:Evaluate: Damu::usage =
88 "Uncertainty of the calculated value of the anomalous magnetic moment of the muon.";
91 "Mixing matrix of the negatively charged charginos.";
94 "Mixing matrix of the positively charged charginos.";
97 "Mixing matrix of the neutralinos.";
99:Evaluate: USt::usage =
100 "Mixing matrix of the stops.";
102:Evaluate: USb::usage =
103 "Mixing matrix of the sbottoms.";
105:Evaluate: UStau::usage =
106 "Mixing matrix of the staus.";
108:Evaluate: USm::usage =
109 "Mixing matrix of the smuons.";
111:Evaluate: MSveL::usage =
112 "Electron-sneutrino mass.";
114:Evaluate: MSvmL::usage =
115 "Muon-sneutrino mass.";
117:Evaluate: MSvtL::usage =
118 "Tau-sneutrino mass.";
120:Evaluate: MSe::usage =
123:Evaluate: MSm::usage =
126:Evaluate: MStau::usage =
129:Evaluate: MSu::usage =
132:Evaluate: MSd::usage =
135:Evaluate: MSc::usage =
138:Evaluate: MSs::usage =
141:Evaluate: MSt::usage =
144:Evaluate: MSb::usage =
147:Evaluate: MChi::usage =
148 "Neutralino masses.";
150:Evaluate: MCha::usage =
153:Evaluate: MhSM::usage =
154 "Standard Model Higgs boson mass.";
156:Evaluate: Mhh::usage =
157 "CP-even Higgs boson masses.";
159:Evaluate: MAh::usage =
160 "CP-odd Higgs boson mass.";
162:Evaluate: MHp::usage =
163 "charged Higgs boson mass.";
165:Evaluate: Mhh::usage =
166 "CP-even Higgs bosons masses.";
168:Evaluate: TB::usage =
169 "tan(beta) = v2/v1 (= vu/vd in the MSSM).";
171:Evaluate: Mu::usage =
172 "Superpotential mu-parameter.";
174:Evaluate: MassB::usage =
175 "Soft-breaking bino mass parameter.";
177:Evaluate: MassWB::usage =
178 "Soft-breaking wino mass parameter.";
180:Evaluate: MassG::usage =
181 "Soft-breaking gluino mass parameter.";
183:Evaluate: mq2::usage =
184 "Soft-breaking squared left-handed squark mass parameters.";
186:Evaluate: mu2::usage =
187 "Soft-breaking squared right-handed up-type squark mass parameters.";
189:Evaluate: md2::usage =
190 "Soft-breaking squared right-handed down-type squark mass parameters.";
192:Evaluate: ml2::usage =
193 "Soft-breaking squared left-handed slepton mass parameters.";
195:Evaluate: me2::usage =
196 "Soft-breaking squared right-handed down-type slepton mass parameters.";
198:Evaluate: Au::usage =
199 "Soft-breaking trilinear couplings between Higgs bosons and up-type squarks.";
201:Evaluate: Ad::usage =
202 "Soft-breaking trilinear couplings between Higgs bosons and down-type squarks.";
204:Evaluate: Ae::usage =
205 "Soft-breaking trilinear couplings between Higgs bosons and down-type sleptons.";
208 "Renormalization scale.";
210:Evaluate: Yu::usage =
211 "Yukawa couplings of up-type quarks.";
213:Evaluate: Yd::usage =
214 "Yukawa couplings of down-type quarks.";
216:Evaluate: Ye::usage =
217 "Yukawa couplings of down-type leptons.";
219:Evaluate: yukawaType::usage =
220 "Yukawa type of the Two-Higgs Doublet Model.";
222:Evaluate: sinBetaMinusAlpha::usage =
223 "Mixing angle sin(beta - alpha) of the Higgs sector in the Two-Higgs Doublet Model.";
225:Evaluate: lambda::usage =
226 "Lagrangian parameters { lambda_1, ... , lambda_7 } in the Two-Higgs Doublet Model.";
228:Evaluate: lambda6::usage =
229 "Lagrangian parameter lambda_6 in the Two-Higgs Doublet Model.";
231:Evaluate: lambda7::usage =
232 "Lagrangian parameter lambda_7 in the Two-Higgs Doublet Model.";
234:Evaluate: m122::usage =
235 "Lagrangian parameter m_{12}^2 in the Two-Higgs Doublet Model.";
237:Evaluate: zetau::usage =
238 "Alignment parameter zeta_u in the Two-Higgs Doublet Model.";
240:Evaluate: zetad::usage =
241 "Alignment parameter zeta_d in the Two-Higgs Doublet Model.";
243:Evaluate: zetal::usage =
244 "Alignment parameter zeta_l in the Two-Higgs Doublet Model.";
246:Evaluate: Deltau::usage =
247 "Yukawa coupling Delta_u in the general Two-Higgs Doublet Model.";
249:Evaluate: Deltad::usage =
250 "Yukawa coupling Delta_d in the general Two-Higgs Doublet Model.";
252:Evaluate: Deltal::usage =
253 "Yukawa coupling Delta_l in the general Two-Higgs Doublet Model.";
255:Evaluate: Piu::usage =
256 "Yukawa coupling Pi_u in the general Two-Higgs Doublet Model.";
258:Evaluate: Pid::usage =
259 "Yukawa coupling Pi_d in the general Two-Higgs Doublet Model.";
261:Evaluate: Pil::usage =
262 "Yukawa coupling Pi_l in the general Two-Higgs Doublet Model.";
264:Evaluate: GM2CalcSetFlags::usage =
265 "GM2CalcSetFlags sets the configuration flags for GM2Calc. Available flags are: {loopOrder, tanBetaResummation, forceOutput, runningCouplings}. Unset flags are set to their default values, see Options[GM2CalcSetFlags]. Use GM2CalcGetFlags[] to retrieve the flags currently set."
267:Evaluate: GM2CalcGetFlags::usage =
268 "GM2CalcGetFlags returns the current configuration flags for GM2Calc."
270:Evaluate: GM2CalcSetSMParameters::usage =
271 "GM2CalcSetSMParameters sets the Standard Model parameters input parameters. Unset parameters are set to their default values, see Options[GM2CalcSetSMParameters]. Use GM2CalcGetSMParameters[] to retrieve the current values of the Standard Model parameters."
273:Evaluate: GM2CalcGetSMParameters::usage =
274 "GM2CalcGetSMParameters returns the Standard Model parameters."
276:Evaluate: GM2CalcAmuSLHAScheme::usage =
277 "GM2CalcAmuSLHAScheme calculates amu and its uncertainty in the MSSM using the given SLHA parameters (SLHA interface). Unset SLHA parameters are set to zero. See Options[GM2CalcAmuSLHAScheme] for all SLHA parameters and their default values."
279:Evaluate: GM2CalcAmuGM2CalcScheme::usage =
280 "GM2CalcAmuGM2CalcScheme calculates amu and its uncertainty in the MSSM using the given parameters in the GM2Calc-specific renormalization scheme (GM2Calc interface). Unset parameters are set to zero. See Options[GM2CalcAmuGM2CalcScheme] for all input parameters in the GM2Calc scheme and their default values."
282:Evaluate: GM2CalcAmuTHDMGaugeBasis::usage =
283 "GM2CalcAmuTHDMGaugeBasis calculates amu and its uncertainty in the Two-Higgs Doublet Model using the given input parameters in the gauge bassis. Unset input parameters are set to zero. See Options[GM2CalcAmuTHDMGaugeBasis] for all parameters and their default values."
285:Evaluate: GM2CalcAmuTHDMMassBasis::usage =
286 "GM2CalcAmuTHDMMassBasis calculates amu and its uncertainty in the Two-Higgs Doublet Model using the given input parameters. Unset input parameters are set to zero. See Options[GM2CalcAmuTHDMMassBasis] for all parameters and their default values."
288:Evaluate: GM2CalcAmuSLHAScheme::error = "`1`";
289:Evaluate: GM2CalcAmuSLHAScheme::warning = "`1`";
291:Evaluate: GM2CalcAmuGM2CalcScheme::error = "`1`";
292:Evaluate: GM2CalcAmuGM2CalcScheme::warning = "`1`";
294:Evaluate: GM2CalcAmuTHDMMassBasis::error = "`1`";
295:Evaluate: GM2CalcAmuTHDMGaugeBasis::error = "`1`";
297:Evaluate: Begin["`Private`"]
300:Function: GM2CalcSetFlags
301:Pattern: GM2CalcSetFlags[OptionsPattern[]]
303 OptionValue[loopOrder],
304 Boole[OptionValue[tanBetaResummation]],
305 Boole[OptionValue[forceOutput]],
306 Boole[OptionValue[runningCouplings]] }
307:ArgumentTypes: {Integer, Integer, Integer, Integer}
311:Evaluate: Options[GM2CalcSetFlags] = {
313 tanBetaResummation -> True,
314 forceOutput -> False,
315 runningCouplings -> False }
317:Evaluate: GM2CalcSetFlags::wronglooporder = "Unsupported loop order: `1`";
320:Function: GM2CalcGetFlags
321:Pattern: GM2CalcGetFlags[]
328:Function: GM2CalcSetSMParameters
329:Pattern: GM2CalcSetSMParameters[OptionsPattern[]]
331 N @ OptionValue[alpha0],
332 N @ OptionValue[alphaMZ],
333 N @ OptionValue[alphaS],
334 N @ OptionValue[MhSM],
338 N @ OptionValue[mcmc],
339 N @ OptionValue[mu2GeV],
340 N @ OptionValue[mbmb],
341 N @ OptionValue[ms2GeV],
342 N @ OptionValue[md2GeV],
343 N @ OptionValue[Mv1],
344 N @ OptionValue[Mv2],
345 N @ OptionValue[Mv3],
349 Re @ N @ OptionValue[CKM][[1,1]],
350 Re @ N @ OptionValue[CKM][[1,2]],
351 Re @ N @ OptionValue[CKM][[1,3]],
352 Re @ N @ OptionValue[CKM][[2,1]],
353 Re @ N @ OptionValue[CKM][[2,2]],
354 Re @ N @ OptionValue[CKM][[2,3]],
355 Re @ N @ OptionValue[CKM][[3,1]],
356 Re @ N @ OptionValue[CKM][[3,2]],
357 Re @ N @ OptionValue[CKM][[3,3]],
358 Im @ N @ OptionValue[CKM][[1,1]],
359 Im @ N @ OptionValue[CKM][[1,2]],
360 Im @ N @ OptionValue[CKM][[1,3]],
361 Im @ N @ OptionValue[CKM][[2,1]],
362 Im @ N @ OptionValue[CKM][[2,2]],
363 Im @ N @ OptionValue[CKM][[2,3]],
364 Im @ N @ OptionValue[CKM][[3,1]],
365 Im @ N @ OptionValue[CKM][[3,2]],
366 Im @ N @ OptionValue[CKM][[3,3]] }
368 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
369 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
370 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
371 Real, Real, Real, Real, Real, Real }
375:Evaluate: Options[GM2CalcSetSMParameters] = {
376 alpha0 -> 0.00729735,
377 alphaMZ -> 0.0077552,
390 ME -> 0.000510998928,
394 CKM -> {{1,0,0}, {0,1,0}, {0,0,1}} }
397:Function: GM2CalcGetSMParameters
398:Pattern: GM2CalcGetSMParameters[]
405:Function: GM2CalcAmuSLHAScheme
406:Pattern: GM2CalcAmuSLHAScheme[OptionsPattern[]]
408 N @ OptionValue[MSvmL],
409 N @ OptionValue[MSm][[1]],
410 N @ OptionValue[MSm][[2]],
411 N @ OptionValue[MChi][[1]],
412 N @ OptionValue[MChi][[2]],
413 N @ OptionValue[MChi][[3]],
414 N @ OptionValue[MChi][[4]],
415 N @ OptionValue[MCha][[1]],
416 N @ OptionValue[MCha][[2]],
417 N @ OptionValue[MAh],
420 N @ OptionValue[MassB],
421 N @ OptionValue[MassWB],
422 N @ OptionValue[MassG],
423 N @ OptionValue[mq2][[1,1]],
424 N @ OptionValue[mq2][[2,2]],
425 N @ OptionValue[mq2][[3,3]],
426 N @ OptionValue[ml2][[1,1]],
427 N @ OptionValue[ml2][[2,2]],
428 N @ OptionValue[ml2][[3,3]],
429 N @ OptionValue[mu2][[1,1]],
430 N @ OptionValue[mu2][[2,2]],
431 N @ OptionValue[mu2][[3,3]],
432 N @ OptionValue[md2][[1,1]],
433 N @ OptionValue[md2][[2,2]],
434 N @ OptionValue[md2][[3,3]],
435 N @ OptionValue[me2][[1,1]],
436 N @ OptionValue[me2][[2,2]],
437 N @ OptionValue[me2][[3,3]],
438 N @ OptionValue[Au][[3,3]],
439 N @ OptionValue[Ad][[3,3]],
440 N @ OptionValue[Ae][[2,2]],
441 N @ OptionValue[Ae][[3,3]],
444 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
445 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
446 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
447 Real, Real, Real, Real, Real }
451:Evaluate: Options[GM2CalcAmuSLHAScheme] = {
454 MChi -> {0, 0, 0, 0},
462 mq2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
463 ml2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
464 mu2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
465 md2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
466 me2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
467 Au -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
468 Ad -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
469 Ae -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
473:Function: GM2CalcAmuGM2CalcScheme
474:Pattern: GM2CalcAmuGM2CalcScheme[OptionsPattern[]]
476 N @ OptionValue[MAh],
479 N @ OptionValue[MassB],
480 N @ OptionValue[MassWB],
481 N @ OptionValue[MassG],
482 N @ OptionValue[mq2][[1,1]],
483 N @ OptionValue[mq2][[2,2]],
484 N @ OptionValue[mq2][[3,3]],
485 N @ OptionValue[ml2][[1,1]],
486 N @ OptionValue[ml2][[2,2]],
487 N @ OptionValue[ml2][[3,3]],
488 N @ OptionValue[mu2][[1,1]],
489 N @ OptionValue[mu2][[2,2]],
490 N @ OptionValue[mu2][[3,3]],
491 N @ OptionValue[md2][[1,1]],
492 N @ OptionValue[md2][[2,2]],
493 N @ OptionValue[md2][[3,3]],
494 N @ OptionValue[me2][[1,1]],
495 N @ OptionValue[me2][[2,2]],
496 N @ OptionValue[me2][[3,3]],
497 N @ OptionValue[Au][[3,3]],
498 N @ OptionValue[Ad][[3,3]],
499 N @ OptionValue[Ae][[2,2]],
500 N @ OptionValue[Ae][[3,3]],
503 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
504 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
505 Real, Real, Real, Real, Real, Real }
509:Evaluate: Options[GM2CalcAmuGM2CalcScheme] = {
516 mq2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
517 ml2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
518 mu2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
519 md2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
520 me2 -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
521 Au -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
522 Ad -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
523 Ae -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
527:Function: GM2CalcAmuTHDMGaugeBasis
528:Pattern: GM2CalcAmuTHDMGaugeBasis[OptionsPattern[]]
530 N @ OptionValue[yukawaType],
531 N @ OptionValue[lambda][[1]],
532 N @ OptionValue[lambda][[2]],
533 N @ OptionValue[lambda][[3]],
534 N @ OptionValue[lambda][[4]],
535 N @ OptionValue[lambda][[5]],
536 N @ OptionValue[lambda][[6]],
537 N @ OptionValue[lambda][[7]],
539 N @ OptionValue[m122],
540 N @ OptionValue[zetau],
541 N @ OptionValue[zetad],
542 N @ OptionValue[zetal],
543 Re @ N @ OptionValue[Deltau][[1,1]],
544 Re @ N @ OptionValue[Deltau][[1,2]],
545 Re @ N @ OptionValue[Deltau][[1,3]],
546 Re @ N @ OptionValue[Deltau][[2,1]],
547 Re @ N @ OptionValue[Deltau][[2,2]],
548 Re @ N @ OptionValue[Deltau][[2,3]],
549 Re @ N @ OptionValue[Deltau][[3,1]],
550 Re @ N @ OptionValue[Deltau][[3,2]],
551 Re @ N @ OptionValue[Deltau][[3,3]],
552 Re @ N @ OptionValue[Deltad][[1,1]],
553 Re @ N @ OptionValue[Deltad][[1,2]],
554 Re @ N @ OptionValue[Deltad][[1,3]],
555 Re @ N @ OptionValue[Deltad][[2,1]],
556 Re @ N @ OptionValue[Deltad][[2,2]],
557 Re @ N @ OptionValue[Deltad][[2,3]],
558 Re @ N @ OptionValue[Deltad][[3,1]],
559 Re @ N @ OptionValue[Deltad][[3,2]],
560 Re @ N @ OptionValue[Deltad][[3,3]],
561 Re @ N @ OptionValue[Deltal][[1,1]],
562 Re @ N @ OptionValue[Deltal][[1,2]],
563 Re @ N @ OptionValue[Deltal][[1,3]],
564 Re @ N @ OptionValue[Deltal][[2,1]],
565 Re @ N @ OptionValue[Deltal][[2,2]],
566 Re @ N @ OptionValue[Deltal][[2,3]],
567 Re @ N @ OptionValue[Deltal][[3,1]],
568 Re @ N @ OptionValue[Deltal][[3,2]],
569 Re @ N @ OptionValue[Deltal][[3,3]],
570 Re @ N @ OptionValue[Piu][[1,1]],
571 Re @ N @ OptionValue[Piu][[1,2]],
572 Re @ N @ OptionValue[Piu][[1,3]],
573 Re @ N @ OptionValue[Piu][[2,1]],
574 Re @ N @ OptionValue[Piu][[2,2]],
575 Re @ N @ OptionValue[Piu][[2,3]],
576 Re @ N @ OptionValue[Piu][[3,1]],
577 Re @ N @ OptionValue[Piu][[3,2]],
578 Re @ N @ OptionValue[Piu][[3,3]],
579 Re @ N @ OptionValue[Pid][[1,1]],
580 Re @ N @ OptionValue[Pid][[1,2]],
581 Re @ N @ OptionValue[Pid][[1,3]],
582 Re @ N @ OptionValue[Pid][[2,1]],
583 Re @ N @ OptionValue[Pid][[2,2]],
584 Re @ N @ OptionValue[Pid][[2,3]],
585 Re @ N @ OptionValue[Pid][[3,1]],
586 Re @ N @ OptionValue[Pid][[3,2]],
587 Re @ N @ OptionValue[Pid][[3,3]],
588 Re @ N @ OptionValue[Pil][[1,1]],
589 Re @ N @ OptionValue[Pil][[1,2]],
590 Re @ N @ OptionValue[Pil][[1,3]],
591 Re @ N @ OptionValue[Pil][[2,1]],
592 Re @ N @ OptionValue[Pil][[2,2]],
593 Re @ N @ OptionValue[Pil][[2,3]],
594 Re @ N @ OptionValue[Pil][[3,1]],
595 Re @ N @ OptionValue[Pil][[3,2]],
596 Re @ N @ OptionValue[Pil][[3,3]] }
598 Integer, Real, Real, Real, Real, Real, Real, Real, Real, Real,
599 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
600 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
601 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
602 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
603 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
604 Real, Real, Real, Real, Real, Real, Real }
608:Evaluate: Options[GM2CalcAmuTHDMGaugeBasis] = {
610 lambda -> { 0, 0, 0, 0, 0, 0, 0 },
616 Deltau -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
617 Deltad -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
618 Deltal -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
619 Piu -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
620 Pid -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
621 Pil -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}} }
624:Function: GM2CalcAmuTHDMMassBasis
625:Pattern: GM2CalcAmuTHDMMassBasis[OptionsPattern[]]
627 N @ OptionValue[yukawaType],
628 N @ OptionValue[Mhh][[1]],
629 N @ OptionValue[Mhh][[2]],
630 N @ OptionValue[MAh],
631 N @ OptionValue[MHp],
632 N @ OptionValue[sinBetaMinusAlpha],
633 N @ OptionValue[lambda6],
634 N @ OptionValue[lambda7],
636 N @ OptionValue[m122],
637 N @ OptionValue[zetau],
638 N @ OptionValue[zetad],
639 N @ OptionValue[zetal],
640 Re @ N @ OptionValue[Deltau][[1,1]],
641 Re @ N @ OptionValue[Deltau][[1,2]],
642 Re @ N @ OptionValue[Deltau][[1,3]],
643 Re @ N @ OptionValue[Deltau][[2,1]],
644 Re @ N @ OptionValue[Deltau][[2,2]],
645 Re @ N @ OptionValue[Deltau][[2,3]],
646 Re @ N @ OptionValue[Deltau][[3,1]],
647 Re @ N @ OptionValue[Deltau][[3,2]],
648 Re @ N @ OptionValue[Deltau][[3,3]],
649 Re @ N @ OptionValue[Deltad][[1,1]],
650 Re @ N @ OptionValue[Deltad][[1,2]],
651 Re @ N @ OptionValue[Deltad][[1,3]],
652 Re @ N @ OptionValue[Deltad][[2,1]],
653 Re @ N @ OptionValue[Deltad][[2,2]],
654 Re @ N @ OptionValue[Deltad][[2,3]],
655 Re @ N @ OptionValue[Deltad][[3,1]],
656 Re @ N @ OptionValue[Deltad][[3,2]],
657 Re @ N @ OptionValue[Deltad][[3,3]],
658 Re @ N @ OptionValue[Deltal][[1,1]],
659 Re @ N @ OptionValue[Deltal][[1,2]],
660 Re @ N @ OptionValue[Deltal][[1,3]],
661 Re @ N @ OptionValue[Deltal][[2,1]],
662 Re @ N @ OptionValue[Deltal][[2,2]],
663 Re @ N @ OptionValue[Deltal][[2,3]],
664 Re @ N @ OptionValue[Deltal][[3,1]],
665 Re @ N @ OptionValue[Deltal][[3,2]],
666 Re @ N @ OptionValue[Deltal][[3,3]],
667 Re @ N @ OptionValue[Piu][[1,1]],
668 Re @ N @ OptionValue[Piu][[1,2]],
669 Re @ N @ OptionValue[Piu][[1,3]],
670 Re @ N @ OptionValue[Piu][[2,1]],
671 Re @ N @ OptionValue[Piu][[2,2]],
672 Re @ N @ OptionValue[Piu][[2,3]],
673 Re @ N @ OptionValue[Piu][[3,1]],
674 Re @ N @ OptionValue[Piu][[3,2]],
675 Re @ N @ OptionValue[Piu][[3,3]],
676 Re @ N @ OptionValue[Pid][[1,1]],
677 Re @ N @ OptionValue[Pid][[1,2]],
678 Re @ N @ OptionValue[Pid][[1,3]],
679 Re @ N @ OptionValue[Pid][[2,1]],
680 Re @ N @ OptionValue[Pid][[2,2]],
681 Re @ N @ OptionValue[Pid][[2,3]],
682 Re @ N @ OptionValue[Pid][[3,1]],
683 Re @ N @ OptionValue[Pid][[3,2]],
684 Re @ N @ OptionValue[Pid][[3,3]],
685 Re @ N @ OptionValue[Pil][[1,1]],
686 Re @ N @ OptionValue[Pil][[1,2]],
687 Re @ N @ OptionValue[Pil][[1,3]],
688 Re @ N @ OptionValue[Pil][[2,1]],
689 Re @ N @ OptionValue[Pil][[2,2]],
690 Re @ N @ OptionValue[Pil][[2,3]],
691 Re @ N @ OptionValue[Pil][[3,1]],
692 Re @ N @ OptionValue[Pil][[3,2]],
693 Re @ N @ OptionValue[Pil][[3,3]] }
695 Integer, Real, Real, Real, Real, Real, Real, Real, Real, Real,
696 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
697 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
698 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
699 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
700 Real, Real, Real, Real, Real, Real, Real, Real, Real, Real,
701 Real, Real, Real, Real, Real, Real, Real }
705:Evaluate: Options[GM2CalcAmuTHDMMassBasis] = {
710 sinBetaMinusAlpha -> 0,
718 Deltau -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
719 Deltad -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
720 Deltal -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
721 Piu -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
722 Pid -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}},
723 Pil -> {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}} }
727:Evaluate: EndPackage[]
730#include "gm2calc/gm2_1loop.h"
731#include "gm2calc/gm2_2loop.h"
732#include "gm2calc/gm2_uncertainty.h"
733#include "gm2calc/gm2_version.h"
734#include "gm2calc/MSSMNoFV_onshell.h"
735#include "gm2calc/SM.h"
736#include "gm2calc/THDM.h"
737#include "gm2_uncertainty_helpers.h"
747#define M_PI 3.14159265358979323846
750#define NELEMS(x) (sizeof(x) / sizeof((x)[0]))
752#define MLPutRule(link,s) \
754 MLPutFunction(link, "Rule", 2); \
755 MLPutSymbol(link, s); \
758#define MLPutRuleToReal(link,v,s) \
760 MLPutRule(link, s); \
761 MLPutReal(link, v); \
764#define MLPutRuleToInteger(link,v,s) \
766 MLPutRule(link, s); \
767 MLPutInteger(link, v); \
770#define MLPutRuleToString(link,v,s) \
772 MLPutRule(link, s); \
773 MLPutString(link, v); \
776#define MLPutRealMatrix(link,v,dim1,dim2) \
778 long dims[] = { dim1, dim2 }; \
779 MLPutDoubleArray(link, v, dims, NULL, NELEMS(dims)); \
782/* macros which access matrix/vector elements through GM2Calc interface
785#define MLPutRealVectorInterface(link,M,dim) \
788 for (unsigned i = 0; i < dim; i++) { \
789 M[i] = gm2calc_mssmnofv_get_##M(model, i); \
791 MLPutRealList(link, M, dim); \
794#define MLPutRealMatrixInterface(link,M,dim1,dim2) \
796 double M[dim1][dim2]; \
797 for (unsigned i = 0; i < dim1; i++) { \
798 for (unsigned k = 0; k < dim2; k++) { \
799 M[i][k] = gm2calc_mssmnofv_get_##M(model, i, k); \
802 MLPutRealMatrix(link, (double*)M, dim1, dim2); \
805#define MLPutComplexMatrixInterface(link,M,dim1,dim2) \
807 MLPutFunction(link, "List", dim1); \
808 for (unsigned i = 0; i < dim1; i++) { \
809 MLPutFunction(link, "List", dim2); \
810 for (unsigned k = 0; k < dim2; k++) { \
811 double re = 0., im = 0.; \
812 re = gm2calc_mssmnofv_get_##M(model, i, k, &im); \
813 MLPutComplex(link, re, im); \
818#define MLPutComplexMatrixTHDM(link,M,dim1,dim2) \
820 MLPutFunction(link, "List", dim1); \
821 for (unsigned i = 0; i < dim1; i++) { \
822 MLPutFunction(link, "List", dim2); \
823 for (unsigned k = 0; k < dim2; k++) { \
824 double re = M##_real[i][k]; \
825 double im = M##_imag[i][k]; \
826 MLPutComplex(link, re, im); \
831#define MLPutRuleToRealVectorInterface(link,M,name,dim) \
833 MLPutFunction(link, "Rule", 2); \
834 MLPutSymbol(link, name); \
835 MLPutRealVectorInterface(link,M,dim); \
838#define MLPutRuleToRealMatrixInterface(link,v,name,dim1,dim2) \
840 MLPutFunction(link, "Rule", 2); \
841 MLPutSymbol(link, name); \
842 MLPutRealMatrixInterface(link,v,dim1,dim2); \
845#define MLPutRuleToComplexMatrixInterface(link,M,name,dim1,dim2) \
847 MLPutFunction(link, "Rule", 2); \
848 MLPutSymbol(link, name); \
849 MLPutComplexMatrixInterface(link,M,dim1,dim2); \
852#define MLPutRuleToComplexMatrixTHDM(link,M,name,dim1,dim2) \
854 MLPutFunction(link, "Rule", 2); \
855 MLPutSymbol(link, name); \
856 MLPutComplexMatrixTHDM(link,M,dim1,dim2); \
859/* global configuration flags */
862 int tanBetaResummation;
864 int runningCouplings;
867 .tanBetaResummation = 1,
869 .runningCouplings = 1
872/* Standard Model parameters */
875/* SUSY parameters for SLHA interface */
876struct SLHA_parameters {
898/* SUSY parameters for GM2Calc interface */
899struct GM2Calc_parameters {
917/******************************************************************/
919void initialize_slha_parameters(struct SLHA_parameters* pars)
922 memset(pars->MSm, 0, sizeof(pars->MSm[0]) * 2);
923 memset(pars->MChi, 0, sizeof(pars->MChi[0]) * 4);
924 memset(pars->MCha, 0, sizeof(pars->MCha[0]) * 2);
931 memset(pars->mq2, 0, sizeof(pars->mq2[0][0]) * 3 * 3);
932 memset(pars->ml2, 0, sizeof(pars->ml2[0][0]) * 3 * 3);
933 memset(pars->mu2, 0, sizeof(pars->mu2[0][0]) * 3 * 3);
934 memset(pars->md2, 0, sizeof(pars->md2[0][0]) * 3 * 3);
935 memset(pars->me2, 0, sizeof(pars->me2[0][0]) * 3 * 3);
936 memset(pars->Au, 0, sizeof(pars->Au[0][0]) * 3 * 3);
937 memset(pars->Ad, 0, sizeof(pars->Ad[0][0]) * 3 * 3);
938 memset(pars->Ae, 0, sizeof(pars->Ae[0][0]) * 3 * 3);
942/******************************************************************/
944void initialize_gm2calc_parameters(struct GM2Calc_parameters* pars)
952 memset(pars->mq2, 0, sizeof(pars->mq2[0][0]) * 3 * 3);
953 memset(pars->ml2, 0, sizeof(pars->ml2[0][0]) * 3 * 3);
954 memset(pars->mu2, 0, sizeof(pars->mu2[0][0]) * 3 * 3);
955 memset(pars->md2, 0, sizeof(pars->md2[0][0]) * 3 * 3);
956 memset(pars->me2, 0, sizeof(pars->me2[0][0]) * 3 * 3);
957 memset(pars->Au, 0, sizeof(pars->Au[0][0]) * 3 * 3);
958 memset(pars->Ad, 0, sizeof(pars->Ad[0][0]) * 3 * 3);
959 memset(pars->Ae, 0, sizeof(pars->Ae[0][0]) * 3 * 3);
963/******************************************************************/
965double sqr(double x) { return x*x; }
967/******************************************************************/
969void MLPutComplex(MLINK link, double re, double im)
974 MLPutFunction(link, "Complex", 2);
980/******************************************************************/
982void put_error_message(const char* function_name,
983 const char* message_tag,
984 const char* message_str)
986 MLPutFunction(stdlink, "CompoundExpression", 2);
987 MLPutFunction(stdlink, "Message", 2);
988 MLPutFunction(stdlink, "MessageName", 2);
989 MLPutSymbol(stdlink, function_name);
990 MLPutString(stdlink, message_tag);
991 MLPutString(stdlink, message_str);
994/******************************************************************/
996gm2calc_error setup_model_slha_scheme(MSSMNoFV_onshell* model,
997 const struct SLHA_parameters* pars)
999 /* fill SM parameters */
1000 gm2calc_mssmnofv_set_alpha_MZ(model, sm.alpha_em_mz);
1001 gm2calc_mssmnofv_set_alpha_thompson(model, sm.alpha_em_0);
1002 gm2calc_mssmnofv_set_g3(model, sqrt(4 * M_PI * sm.alpha_s_mz));
1003 gm2calc_mssmnofv_set_MT_pole(model, sm.mu[2]);
1004 gm2calc_mssmnofv_set_MB_running(model, sm.md[2]);
1005 gm2calc_mssmnofv_set_MM_pole(model, sm.ml[1]);
1006 gm2calc_mssmnofv_set_ML_pole(model, sm.ml[2]);
1007 gm2calc_mssmnofv_set_MW_pole(model, sm.mw);
1008 gm2calc_mssmnofv_set_MZ_pole(model, sm.mz);
1010 /* fill pole masses */
1011 gm2calc_mssmnofv_set_MSvmL_pole(model, pars->MSvmL);
1012 gm2calc_mssmnofv_set_MSm_pole(model, 0, pars->MSm[0]);
1013 gm2calc_mssmnofv_set_MSm_pole(model, 1, pars->MSm[1]);
1014 gm2calc_mssmnofv_set_MChi_pole(model, 0, pars->MChi[0]);
1015 gm2calc_mssmnofv_set_MChi_pole(model, 1, pars->MChi[1]);
1016 gm2calc_mssmnofv_set_MChi_pole(model, 2, pars->MChi[2]);
1017 gm2calc_mssmnofv_set_MChi_pole(model, 3, pars->MChi[3]);
1018 gm2calc_mssmnofv_set_MCha_pole(model, 0, pars->MCha[0]);
1019 gm2calc_mssmnofv_set_MCha_pole(model, 1, pars->MCha[1]);
1020 gm2calc_mssmnofv_set_MAh_pole(model, pars->MAh);
1022 /* fill DR-bar parameters */
1023 gm2calc_mssmnofv_set_TB(model, pars->TB);
1024 gm2calc_mssmnofv_set_Mu(model, pars->Mu);
1025 gm2calc_mssmnofv_set_MassB(model, pars->MassB);
1026 gm2calc_mssmnofv_set_MassWB(model, pars->MassWB);
1027 gm2calc_mssmnofv_set_MassG(model, pars->MassG);
1028 for (unsigned i = 0; i < 3; i++) {
1029 for (unsigned k = 0; k < 3; k++) {
1030 gm2calc_mssmnofv_set_mq2(model, i, k, pars->mq2[i][k]);
1031 gm2calc_mssmnofv_set_ml2(model, i, k, pars->ml2[i][k]);
1032 gm2calc_mssmnofv_set_mu2(model, i, k, pars->mu2[i][k]);
1033 gm2calc_mssmnofv_set_md2(model, i, k, pars->md2[i][k]);
1034 gm2calc_mssmnofv_set_me2(model, i, k, pars->me2[i][k]);
1035 gm2calc_mssmnofv_set_Au(model, i, k, pars->Au[i][k]);
1036 gm2calc_mssmnofv_set_Ad(model, i, k, pars->Ad[i][k]);
1037 gm2calc_mssmnofv_set_Ae(model, i, k, pars->Ae[i][k]);
1040 gm2calc_mssmnofv_set_scale(model, pars->Q);
1042 /* convert DR-bar parameters to on-shell */
1043 const gm2calc_error error = gm2calc_mssmnofv_convert_to_onshell(model);
1048/******************************************************************/
1050gm2calc_error setup_model_gm2calc_scheme(MSSMNoFV_onshell* model,
1051 const struct GM2Calc_parameters* pars)
1053 /* fill SM parameters */
1054 gm2calc_mssmnofv_set_alpha_MZ(model, sm.alpha_em_mz);
1055 gm2calc_mssmnofv_set_alpha_thompson(model, sm.alpha_em_0);
1056 gm2calc_mssmnofv_set_g3(model, sqrt(4 * M_PI * sm.alpha_s_mz));
1057 gm2calc_mssmnofv_set_MT_pole(model, sm.mu[2]);
1058 gm2calc_mssmnofv_set_MB_running(model, sm.md[2]);
1059 gm2calc_mssmnofv_set_MM_pole(model, sm.ml[1]);
1060 gm2calc_mssmnofv_set_ML_pole(model, sm.ml[2]);
1061 gm2calc_mssmnofv_set_MW_pole(model, sm.mw);
1062 gm2calc_mssmnofv_set_MZ_pole(model, sm.mz);
1064 /* fill DR-bar parameters */
1065 gm2calc_mssmnofv_set_TB(model, pars->TB);
1067 /* fill on-shell parameters */
1068 gm2calc_mssmnofv_set_Mu(model, pars->Mu);
1069 gm2calc_mssmnofv_set_MassB(model, pars->MassB);
1070 gm2calc_mssmnofv_set_MassWB(model, pars->MassWB);
1071 gm2calc_mssmnofv_set_MassG(model, pars->MassG);
1072 gm2calc_mssmnofv_set_MAh_pole(model, pars->MAh);
1073 gm2calc_mssmnofv_set_scale(model, pars->Q);
1075 /* fill DR-bar parameters */
1076 for (unsigned i = 0; i < 3; i++) {
1077 for (unsigned k = 0; k < 3; k++) {
1078 gm2calc_mssmnofv_set_mq2(model, i, k, pars->mq2[i][k]);
1079 gm2calc_mssmnofv_set_ml2(model, i, k, pars->ml2[i][k]);
1080 gm2calc_mssmnofv_set_mu2(model, i, k, pars->mu2[i][k]);
1081 gm2calc_mssmnofv_set_md2(model, i, k, pars->md2[i][k]);
1082 gm2calc_mssmnofv_set_me2(model, i, k, pars->me2[i][k]);
1083 gm2calc_mssmnofv_set_Au(model, i, k, pars->Au[i][k]);
1084 gm2calc_mssmnofv_set_Ad(model, i, k, pars->Ad[i][k]);
1085 gm2calc_mssmnofv_set_Ae(model, i, k, pars->Ae[i][k]);
1089 /* convert DR-bar parameters to on-shell */
1090 const gm2calc_error error = gm2calc_mssmnofv_calculate_masses(model);
1095/******************************************************************/
1097void calculate_amu_mssmnofv(MSSMNoFV_onshell* model, double* amu1L, double* amu2L)
1099 if (config_flags.loopOrder > 0) {
1100 if (config_flags.tanBetaResummation) {
1101 *amu1L = gm2calc_mssmnofv_calculate_amu_1loop(model);
1103 *amu1L = gm2calc_mssmnofv_calculate_amu_1loop_non_tan_beta_resummed(model);
1106 if (config_flags.loopOrder > 1) {
1107 if (config_flags.tanBetaResummation) {
1108 *amu2L = gm2calc_mssmnofv_calculate_amu_2loop(model);
1110 *amu2L = gm2calc_mssmnofv_calculate_amu_2loop_non_tan_beta_resummed(model);
1115/******************************************************************/
1117double calculate_uncertainty_mssmnofv(MSSMNoFV_onshell* model, double amu1L, double amu2L)
1119 double delta_amu = 0.;
1121 if (config_flags.loopOrder == 0) {
1122 delta_amu = gm2calc_mssmnofv_calculate_uncertainty_amu_0loop_amu1L(model, amu1L);
1123 } else if (config_flags.loopOrder == 1) {
1124 delta_amu = gm2calc_mssmnofv_calculate_uncertainty_amu_1loop_amu2L(model, amu2L);
1125 } else if (config_flags.loopOrder > 1) {
1126 delta_amu = gm2calc_mssmnofv_calculate_uncertainty_amu_2loop(model);
1132/******************************************************************/
1134void calculate_amu_thdm(gm2calc_THDM* model, double* amu1L, double* amu2LF, double* amu2LB)
1136 if (config_flags.loopOrder > 0) {
1137 *amu1L = gm2calc_thdm_calculate_amu_1loop(model);
1139 if (config_flags.loopOrder > 1) {
1140 *amu2LF = gm2calc_thdm_calculate_amu_2loop_fermionic(model);
1141 *amu2LB = gm2calc_thdm_calculate_amu_2loop_bosonic(model);
1145/******************************************************************/
1147double calculate_uncertainty_thdm(gm2calc_THDM* model, double amu1L, double amu2L)
1149 double delta_amu = 0.;
1151 if (config_flags.loopOrder == 0) {
1152 delta_amu = gm2calc_thdm_calculate_uncertainty_amu_0loop_amu1L_amu2L(model, amu1L, amu2L);
1153 } else if (config_flags.loopOrder == 1) {
1154 delta_amu = gm2calc_thdm_calculate_uncertainty_amu_1loop_amu1L_amu2L(model, amu1L, amu2L);
1155 } else if (config_flags.loopOrder > 1) {
1156 delta_amu = gm2calc_thdm_calculate_uncertainty_amu_2loop_amu1L_amu2L(model, amu1L, amu2L);
1162/******************************************************************/
1164static void print_package()
1166 static int do_print = 1;
1169 printf("===========================\n");
1170 printf("GM2Calc " GM2CALC_VERSION "\n");
1171 printf("http://gm2calc.hepforge.org\n");
1172 printf("===========================\n");
1178/******************************************************************/
1180int GM2CalcSetFlags(int loopOrder_, int tanBetaResummation_, int forceOutput_,
1181 int runningCouplings_)
1183 char loop_order_str[12];
1187 if (loopOrder_ < 0 || loopOrder_ > 2) {
1188 snprintf(loop_order_str, sizeof(loop_order_str), "%d", loopOrder_);
1189 put_error_message("GM2CalcSetFlags", "wronglooporder", loop_order_str);
1192 config_flags.loopOrder = loopOrder_;
1193 config_flags.tanBetaResummation = tanBetaResummation_;
1194 config_flags.forceOutput = forceOutput_;
1195 config_flags.runningCouplings = runningCouplings_;
1200/******************************************************************/
1202void GM2CalcGetFlags(void)
1204 MLPutFunction(stdlink, "List", 4);
1205 MLPutRuleToInteger(stdlink, config_flags.loopOrder, "loopOrder");
1206 MLPutRuleToString(stdlink, config_flags.tanBetaResummation ? "True" : "False", "tanBetaResummation");
1207 MLPutRuleToString(stdlink, config_flags.forceOutput ? "True" : "False", "forceOutput");
1208 MLPutRuleToString(stdlink, config_flags.runningCouplings ? "True" : "False", "runningCouplings");
1209 MLEndPacket(stdlink);
1212/******************************************************************/
1214int GM2CalcSetSMParameters(
1233 double CKM_real_11_,
1234 double CKM_real_12_,
1235 double CKM_real_13_,
1236 double CKM_real_21_,
1237 double CKM_real_22_,
1238 double CKM_real_23_,
1239 double CKM_real_31_,
1240 double CKM_real_32_,
1241 double CKM_real_33_,
1242 double CKM_imag_11_,
1243 double CKM_imag_12_,
1244 double CKM_imag_13_,
1245 double CKM_imag_21_,
1246 double CKM_imag_22_,
1247 double CKM_imag_23_,
1248 double CKM_imag_31_,
1249 double CKM_imag_32_,
1250 double CKM_imag_33_)
1252 sm.alpha_em_0 = alpha0_;
1253 sm.alpha_em_mz = alphaMZ_;
1254 sm.alpha_s_mz = alphaS_;
1270 sm.ckm_real[0][0] = CKM_real_11_;
1271 sm.ckm_real[0][1] = CKM_real_12_;
1272 sm.ckm_real[0][2] = CKM_real_13_;
1273 sm.ckm_real[1][0] = CKM_real_21_;
1274 sm.ckm_real[1][1] = CKM_real_22_;
1275 sm.ckm_real[1][2] = CKM_real_23_;
1276 sm.ckm_real[2][0] = CKM_real_31_;
1277 sm.ckm_real[2][1] = CKM_real_32_;
1278 sm.ckm_real[2][2] = CKM_real_33_;
1279 sm.ckm_imag[0][0] = CKM_imag_11_;
1280 sm.ckm_imag[0][1] = CKM_imag_12_;
1281 sm.ckm_imag[0][2] = CKM_imag_13_;
1282 sm.ckm_imag[1][0] = CKM_imag_21_;
1283 sm.ckm_imag[1][1] = CKM_imag_22_;
1284 sm.ckm_imag[1][2] = CKM_imag_23_;
1285 sm.ckm_imag[2][0] = CKM_imag_31_;
1286 sm.ckm_imag[2][1] = CKM_imag_32_;
1287 sm.ckm_imag[2][2] = CKM_imag_33_;
1292/******************************************************************/
1294void GM2CalcGetSMParameters(void)
1296 MLPutFunction(stdlink, "List", 19);
1298 MLPutRuleToReal(stdlink, sm.alpha_em_0, "alpha0");
1299 MLPutRuleToReal(stdlink, sm.alpha_em_mz, "alphaMZ");
1300 MLPutRuleToReal(stdlink, sm.alpha_s_mz, "alphaS");
1301 MLPutRuleToReal(stdlink, sm.mh, "MhSM");
1302 MLPutRuleToReal(stdlink, sm.mw, "MW");
1303 MLPutRuleToReal(stdlink, sm.mz, "MZ");
1304 MLPutRuleToReal(stdlink, sm.mu[2], "MT");
1305 MLPutRuleToReal(stdlink, sm.mu[1], "mcmc");
1306 MLPutRuleToReal(stdlink, sm.mu[0], "mu2GeV");
1307 MLPutRuleToReal(stdlink, sm.md[2], "mbmb");
1308 MLPutRuleToReal(stdlink, sm.md[1], "ms2GeV");
1309 MLPutRuleToReal(stdlink, sm.md[0], "md2GeV");
1310 MLPutRuleToReal(stdlink, sm.mv[2], "Mv1");
1311 MLPutRuleToReal(stdlink, sm.mv[1], "Mv2");
1312 MLPutRuleToReal(stdlink, sm.mv[0], "Mv3");
1313 MLPutRuleToReal(stdlink, sm.ml[2], "ML");
1314 MLPutRuleToReal(stdlink, sm.ml[1], "MM");
1315 MLPutRuleToReal(stdlink, sm.ml[0], "ME");
1316 MLPutRuleToComplexMatrixTHDM(stdlink, sm.ckm, "CKM", 3, 3);
1318 MLEndPacket(stdlink);
1321/******************************************************************/
1323void create_result_list(MSSMNoFV_onshell* model)
1325 double amu1L = 0, amu2L = 0;
1326 calculate_amu_mssmnofv(model, &amu1L, &amu2L);
1327 const double amu = amu1L + amu2L;
1328 const double Damu = calculate_uncertainty_mssmnofv(model, amu1L, amu2L);
1330 MLPutFunction(stdlink, "List", 47);
1332 MLPutRuleToReal(stdlink, amu, "amu");
1333 MLPutRuleToReal(stdlink, amu1L, "amu1L");
1334 MLPutRuleToReal(stdlink, amu2L, "amu2L");
1335 MLPutRuleToReal(stdlink, Damu, "Damu");
1337 MLPutRuleToReal(stdlink, sqr(gm2calc_mssmnofv_get_EL(model))/(4.*M_PI), "alphaMZ");
1338 MLPutRuleToReal(stdlink, sqr(gm2calc_mssmnofv_get_EL0(model))/(4.*M_PI), "alpha0");
1339 MLPutRuleToReal(stdlink, sqr(gm2calc_mssmnofv_get_g3(model))/(4.*M_PI), "alphaS");
1340 /* on-shell masses and parameters [40] */
1341 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MM(model), "MM");
1342 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MT(model), "MT");
1343 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MBMB(model), "mbmb");
1344 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MB(model), "mbMZ");
1345 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_ML(model), "ML");
1346 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MW(model), "MW");
1347 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MZ(model), "MZ");
1348 MLPutRuleToRealVectorInterface(stdlink, MSm, "MSm", 2);
1349 MLPutRuleToRealMatrixInterface(stdlink, USm, "USm", 2, 2);
1350 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MSvmL(model), "MSvmL");
1351 MLPutRuleToRealVectorInterface(stdlink, MStau, "MStau", 2);
1352 MLPutRuleToRealMatrixInterface(stdlink, UStau, "UStau", 2, 2);
1353 MLPutRuleToRealVectorInterface(stdlink, MSt, "MSt", 2);
1354 MLPutRuleToRealMatrixInterface(stdlink, USt, "USt", 2, 2);
1355 MLPutRuleToRealVectorInterface(stdlink, MSb, "MSb", 2);
1356 MLPutRuleToRealMatrixInterface(stdlink, USb, "USb", 2, 2);
1357 MLPutRuleToRealVectorInterface(stdlink, MCha, "MCha", 2);
1358 MLPutRuleToComplexMatrixInterface(stdlink, UM, "UM", 2, 2);
1359 MLPutRuleToComplexMatrixInterface(stdlink, UP, "UP", 2, 2);
1360 MLPutRuleToRealVectorInterface(stdlink, MChi, "MChi", 4);
1361 MLPutRuleToComplexMatrixInterface(stdlink, ZN, "ZN", 4, 4);
1362 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MAh(model), "MAh");
1363 MLPutRuleToRealVectorInterface(stdlink, Mhh, "Mhh", 2);
1364 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_TB(model), "TB");
1365 MLPutRuleToRealMatrixInterface(stdlink, Yu, "Yu", 3, 3);
1366 MLPutRuleToRealMatrixInterface(stdlink, Yd, "Yd", 3, 3);
1367 MLPutRuleToRealMatrixInterface(stdlink, Ye, "Ye", 3, 3);
1368 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_Mu(model), "Mu");
1369 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MassB(model), "MassB");
1370 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MassWB(model), "MassWB");
1371 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_MassG(model), "MassG");
1372 MLPutRuleToRealMatrixInterface(stdlink, mq2, "mq2", 3, 3);
1373 MLPutRuleToRealMatrixInterface(stdlink, ml2, "ml2", 3, 3);
1374 MLPutRuleToRealMatrixInterface(stdlink, mu2, "mu2", 3, 3);
1375 MLPutRuleToRealMatrixInterface(stdlink, md2, "md2", 3, 3);
1376 MLPutRuleToRealMatrixInterface(stdlink, me2, "me2", 3, 3);
1377 MLPutRuleToRealMatrixInterface(stdlink, Au, "Au", 3, 3);
1378 MLPutRuleToRealMatrixInterface(stdlink, Ad, "Ad", 3, 3);
1379 MLPutRuleToRealMatrixInterface(stdlink, Ae, "Ae", 3, 3);
1380 MLPutRuleToReal(stdlink, gm2calc_mssmnofv_get_scale(model), "Q");
1381 MLEndPacket(stdlink);
1384/******************************************************************/
1386void create_error_output(void)
1388 MLPutFunction(stdlink, "List", 0);
1389 MLEndPacket(stdlink);
1392/******************************************************************/
1394void GM2CalcAmuSLHAScheme(
1432 struct SLHA_parameters pars;
1433 initialize_slha_parameters(&pars);
1435 pars.MSvmL = MSvmL_;
1436 pars.MSm[0] = MSm_1_;
1437 pars.MSm[1] = MSm_2_;
1438 pars.MChi[0] = MChi_1_;
1439 pars.MChi[1] = MChi_2_;
1440 pars.MChi[2] = MChi_3_;
1441 pars.MChi[3] = MChi_4_;
1442 pars.MCha[0] = MCha_1_;
1443 pars.MCha[1] = MCha_2_;
1447 pars.MassB = MassB_;
1448 pars.MassWB = MassWB_;
1449 pars.MassG = MassG_;
1450 pars.mq2[0][0] = mq2_11_;
1451 pars.mq2[1][1] = mq2_22_;
1452 pars.mq2[2][2] = mq2_33_;
1453 pars.ml2[0][0] = ml2_11_;
1454 pars.ml2[1][1] = ml2_22_;
1455 pars.ml2[2][2] = ml2_33_;
1456 pars.mu2[0][0] = mu2_11_;
1457 pars.mu2[1][1] = mu2_22_;
1458 pars.mu2[2][2] = mu2_33_;
1459 pars.md2[0][0] = md2_11_;
1460 pars.md2[1][1] = md2_22_;
1461 pars.md2[2][2] = md2_33_;
1462 pars.me2[0][0] = me2_11_;
1463 pars.me2[1][1] = me2_22_;
1464 pars.me2[2][2] = me2_33_;
1465 pars.Au[2][2] = Au_33_;
1466 pars.Ad[2][2] = Ad_33_;
1467 pars.Ae[1][1] = Ae_22_;
1468 pars.Ae[2][2] = Ae_33_;
1471 MSSMNoFV_onshell* model = gm2calc_mssmnofv_new();
1472 const gm2calc_error error = setup_model_slha_scheme(model, &pars);
1474 if (gm2calc_mssmnofv_have_warning(model)) {
1476 gm2calc_mssmnofv_get_warnings(model, msg, sizeof(msg));
1477 put_error_message("GM2CalcAmuSLHAScheme", "warning", msg);
1480 if (error != gm2calc_NoError) {
1481 put_error_message("GM2CalcAmuSLHAScheme", "error",
1482 gm2calc_error_str(error));
1485 if (gm2calc_mssmnofv_have_problem(model)) {
1487 gm2calc_mssmnofv_get_problems(model, msg, sizeof(msg));
1488 put_error_message("GM2CalcAmuSLHAScheme", "error", msg);
1491 if (error == gm2calc_NoError || config_flags.forceOutput) {
1492 create_result_list(model);
1494 create_error_output();
1497 gm2calc_mssmnofv_free(model);
1500/******************************************************************/
1502void GM2CalcAmuGM2CalcScheme(
1530 struct GM2Calc_parameters pars;
1531 initialize_gm2calc_parameters(&pars);
1536 pars.MassB = MassB_;
1537 pars.MassWB = MassWB_;
1538 pars.MassG = MassG_;
1539 pars.mq2[0][0] = mq2_11_;
1540 pars.mq2[1][1] = mq2_22_;
1541 pars.mq2[2][2] = mq2_33_;
1542 pars.ml2[0][0] = ml2_11_;
1543 pars.ml2[1][1] = ml2_22_;
1544 pars.ml2[2][2] = ml2_33_;
1545 pars.mu2[0][0] = mu2_11_;
1546 pars.mu2[1][1] = mu2_22_;
1547 pars.mu2[2][2] = mu2_33_;
1548 pars.md2[0][0] = md2_11_;
1549 pars.md2[1][1] = md2_22_;
1550 pars.md2[2][2] = md2_33_;
1551 pars.me2[0][0] = me2_11_;
1552 pars.me2[1][1] = me2_22_;
1553 pars.me2[2][2] = me2_33_;
1554 pars.Au[2][2] = Au_33_;
1555 pars.Ad[2][2] = Ad_33_;
1556 pars.Ae[1][1] = Ae_22_;
1557 pars.Ae[2][2] = Ae_33_;
1560 MSSMNoFV_onshell* model = gm2calc_mssmnofv_new();
1561 const gm2calc_error error = setup_model_gm2calc_scheme(model, &pars);
1563 if (gm2calc_mssmnofv_have_warning(model)) {
1565 gm2calc_mssmnofv_get_warnings(model, msg, sizeof(msg));
1566 put_error_message("GM2CalcAmuGM2CalcScheme", "warning", msg);
1569 if (error != gm2calc_NoError) {
1570 put_error_message("GM2CalcAmuGM2CalcScheme", "error",
1571 gm2calc_error_str(error));
1574 if (gm2calc_mssmnofv_have_problem(model)) {
1576 gm2calc_mssmnofv_get_problems(model, msg, sizeof(msg));
1577 put_error_message("GM2CalcAmuGM2CalcScheme", "error", msg);
1580 if (error == gm2calc_NoError || config_flags.forceOutput) {
1581 create_result_list(model);
1583 create_error_output();
1586 gm2calc_mssmnofv_free(model);
1589/******************************************************************/
1591void GM2CalcAmuTHDMGaugeBasis(
1660 if (yukawa_type_ < 1 || yukawa_type_ > 6) {
1661 put_error_message("GM2CalcAmuTHDMGaugeBasis", "error",
1662 "yukawaType must be between 1 and 6.");
1663 create_error_output();
1667 gm2calc_THDM_gauge_basis basis;
1668 basis.yukawa_type = int_to_c_yukawa_type(yukawa_type_);
1669 basis.lambda[0] = lambda_1_;
1670 basis.lambda[1] = lambda_2_;
1671 basis.lambda[2] = lambda_3_;
1672 basis.lambda[3] = lambda_4_;
1673 basis.lambda[4] = lambda_5_;
1674 basis.lambda[5] = lambda_6_;
1675 basis.lambda[6] = lambda_7_;
1676 basis.tan_beta = TB_;
1678 basis.zeta_u = zeta_u_;
1679 basis.zeta_d = zeta_d_;
1680 basis.zeta_l = zeta_l_;
1681 basis.Delta_u[0][0] = Deltau_11_;
1682 basis.Delta_u[0][1] = Deltau_12_;
1683 basis.Delta_u[0][2] = Deltau_13_;
1684 basis.Delta_u[1][0] = Deltau_21_;
1685 basis.Delta_u[1][1] = Deltau_22_;
1686 basis.Delta_u[1][2] = Deltau_23_;
1687 basis.Delta_u[2][0] = Deltau_31_;
1688 basis.Delta_u[2][1] = Deltau_32_;
1689 basis.Delta_u[2][2] = Deltau_33_;
1690 basis.Delta_d[0][0] = Deltad_11_;
1691 basis.Delta_d[0][1] = Deltad_12_;
1692 basis.Delta_d[0][2] = Deltad_13_;
1693 basis.Delta_d[1][0] = Deltad_21_;
1694 basis.Delta_d[1][1] = Deltad_22_;
1695 basis.Delta_d[1][2] = Deltad_23_;
1696 basis.Delta_d[2][0] = Deltad_31_;
1697 basis.Delta_d[2][1] = Deltad_32_;
1698 basis.Delta_d[2][2] = Deltad_33_;
1699 basis.Delta_l[0][0] = Deltal_11_;
1700 basis.Delta_l[0][1] = Deltal_12_;
1701 basis.Delta_l[0][2] = Deltal_13_;
1702 basis.Delta_l[1][0] = Deltal_21_;
1703 basis.Delta_l[1][1] = Deltal_22_;
1704 basis.Delta_l[1][2] = Deltal_23_;
1705 basis.Delta_l[2][0] = Deltal_31_;
1706 basis.Delta_l[2][1] = Deltal_32_;
1707 basis.Delta_l[2][2] = Deltal_33_;
1708 basis.Pi_u[0][0] = Piu_11_;
1709 basis.Pi_u[0][1] = Piu_12_;
1710 basis.Pi_u[0][2] = Piu_13_;
1711 basis.Pi_u[1][0] = Piu_21_;
1712 basis.Pi_u[1][1] = Piu_22_;
1713 basis.Pi_u[1][2] = Piu_23_;
1714 basis.Pi_u[2][0] = Piu_31_;
1715 basis.Pi_u[2][1] = Piu_32_;
1716 basis.Pi_u[2][2] = Piu_33_;
1717 basis.Pi_d[0][0] = Pid_11_;
1718 basis.Pi_d[0][1] = Pid_12_;
1719 basis.Pi_d[0][2] = Pid_13_;
1720 basis.Pi_d[1][0] = Pid_21_;
1721 basis.Pi_d[1][1] = Pid_22_;
1722 basis.Pi_d[1][2] = Pid_23_;
1723 basis.Pi_d[2][0] = Pid_31_;
1724 basis.Pi_d[2][1] = Pid_32_;
1725 basis.Pi_d[2][2] = Pid_33_;
1726 basis.Pi_l[0][0] = Pil_11_;
1727 basis.Pi_l[0][1] = Pil_12_;
1728 basis.Pi_l[0][2] = Pil_13_;
1729 basis.Pi_l[1][0] = Pil_21_;
1730 basis.Pi_l[1][1] = Pil_22_;
1731 basis.Pi_l[1][2] = Pil_23_;
1732 basis.Pi_l[2][0] = Pil_31_;
1733 basis.Pi_l[2][1] = Pil_32_;
1734 basis.Pi_l[2][2] = Pil_33_;
1736 gm2calc_THDM_config config;
1737 config.force_output = config_flags.forceOutput;
1738 config.running_couplings = config_flags.runningCouplings;
1740 gm2calc_THDM* model = 0;
1741 gm2calc_error error = gm2calc_thdm_new_with_gauge_basis(&model, &basis, &sm, &config);
1743 if (error == gm2calc_NoError) {
1744 double amu1L = 0, amu2LF = 0, amu2LB = 0;
1745 calculate_amu_thdm(model, &amu1L, &amu2LF, &amu2LB);
1746 const double damu = calculate_uncertainty_thdm(model, amu1L, amu2LF + amu2LB);
1748 MLPutFunction(stdlink, "List", 5);
1749 MLPutRuleToReal(stdlink, amu1L + amu2LF + amu2LB, "amu");
1750 MLPutRuleToReal(stdlink, amu1L, "amu1L");
1751 MLPutRuleToReal(stdlink, amu2LF, "amu2LF");
1752 MLPutRuleToReal(stdlink, amu2LB, "amu2LB");
1753 MLPutRuleToReal(stdlink, damu, "Damu");
1754 MLEndPacket(stdlink);
1756 put_error_message("GM2CalcAmuTHDMMassBasis", "error",
1757 gm2calc_error_str(error));
1758 create_error_output();
1761 gm2calc_thdm_free(model);
1764/******************************************************************/
1766void GM2CalcAmuTHDMMassBasis(
1772 double sin_beta_minus_alpha_,
1835 if (yukawa_type_ < 1 || yukawa_type_ > 6) {
1836 put_error_message("GM2CalcAmuTHDMMassBasis", "error",
1837 "yukawaType must be between 1 and 6.");
1838 create_error_output();
1842 gm2calc_THDM_mass_basis basis;
1843 basis.yukawa_type = int_to_c_yukawa_type(yukawa_type_);
1848 basis.sin_beta_minus_alpha = sin_beta_minus_alpha_;
1849 basis.lambda_6 = lambda_6_;
1850 basis.lambda_7 = lambda_7_;
1851 basis.tan_beta = TB_;
1853 basis.zeta_u = zeta_u_;
1854 basis.zeta_d = zeta_d_;
1855 basis.zeta_l = zeta_l_;
1856 basis.Delta_u[0][0] = Deltau_11_;
1857 basis.Delta_u[0][1] = Deltau_12_;
1858 basis.Delta_u[0][2] = Deltau_13_;
1859 basis.Delta_u[1][0] = Deltau_21_;
1860 basis.Delta_u[1][1] = Deltau_22_;
1861 basis.Delta_u[1][2] = Deltau_23_;
1862 basis.Delta_u[2][0] = Deltau_31_;
1863 basis.Delta_u[2][1] = Deltau_32_;
1864 basis.Delta_u[2][2] = Deltau_33_;
1865 basis.Delta_d[0][0] = Deltad_11_;
1866 basis.Delta_d[0][1] = Deltad_12_;
1867 basis.Delta_d[0][2] = Deltad_13_;
1868 basis.Delta_d[1][0] = Deltad_21_;
1869 basis.Delta_d[1][1] = Deltad_22_;
1870 basis.Delta_d[1][2] = Deltad_23_;
1871 basis.Delta_d[2][0] = Deltad_31_;
1872 basis.Delta_d[2][1] = Deltad_32_;
1873 basis.Delta_d[2][2] = Deltad_33_;
1874 basis.Delta_l[0][0] = Deltal_11_;
1875 basis.Delta_l[0][1] = Deltal_12_;
1876 basis.Delta_l[0][2] = Deltal_13_;
1877 basis.Delta_l[1][0] = Deltal_21_;
1878 basis.Delta_l[1][1] = Deltal_22_;
1879 basis.Delta_l[1][2] = Deltal_23_;
1880 basis.Delta_l[2][0] = Deltal_31_;
1881 basis.Delta_l[2][1] = Deltal_32_;
1882 basis.Delta_l[2][2] = Deltal_33_;
1883 basis.Pi_u[0][0] = Piu_11_;
1884 basis.Pi_u[0][1] = Piu_12_;
1885 basis.Pi_u[0][2] = Piu_13_;
1886 basis.Pi_u[1][0] = Piu_21_;
1887 basis.Pi_u[1][1] = Piu_22_;
1888 basis.Pi_u[1][2] = Piu_23_;
1889 basis.Pi_u[2][0] = Piu_31_;
1890 basis.Pi_u[2][1] = Piu_32_;
1891 basis.Pi_u[2][2] = Piu_33_;
1892 basis.Pi_d[0][0] = Pid_11_;
1893 basis.Pi_d[0][1] = Pid_12_;
1894 basis.Pi_d[0][2] = Pid_13_;
1895 basis.Pi_d[1][0] = Pid_21_;
1896 basis.Pi_d[1][1] = Pid_22_;
1897 basis.Pi_d[1][2] = Pid_23_;
1898 basis.Pi_d[2][0] = Pid_31_;
1899 basis.Pi_d[2][1] = Pid_32_;
1900 basis.Pi_d[2][2] = Pid_33_;
1901 basis.Pi_l[0][0] = Pil_11_;
1902 basis.Pi_l[0][1] = Pil_12_;
1903 basis.Pi_l[0][2] = Pil_13_;
1904 basis.Pi_l[1][0] = Pil_21_;
1905 basis.Pi_l[1][1] = Pil_22_;
1906 basis.Pi_l[1][2] = Pil_23_;
1907 basis.Pi_l[2][0] = Pil_31_;
1908 basis.Pi_l[2][1] = Pil_32_;
1909 basis.Pi_l[2][2] = Pil_33_;
1911 gm2calc_THDM_config config;
1912 config.force_output = config_flags.forceOutput;
1913 config.running_couplings = config_flags.runningCouplings;
1915 gm2calc_THDM* model = 0;
1916 gm2calc_error error = gm2calc_thdm_new_with_mass_basis(&model, &basis, &sm, &config);
1918 if (error == gm2calc_NoError) {
1919 double amu1L = 0, amu2LF = 0, amu2LB = 0;
1920 calculate_amu_thdm(model, &amu1L, &amu2LF, &amu2LB);
1921 const double damu = calculate_uncertainty_thdm(model, amu1L, amu2LF + amu2LB);
1923 MLPutFunction(stdlink, "List", 5);
1924 MLPutRuleToReal(stdlink, amu1L + amu2LF + amu2LB, "amu");
1925 MLPutRuleToReal(stdlink, amu1L, "amu1L");
1926 MLPutRuleToReal(stdlink, amu2LF, "amu2LF");
1927 MLPutRuleToReal(stdlink, amu2LB, "amu2LB");
1928 MLPutRuleToReal(stdlink, damu, "Damu");
1929 MLEndPacket(stdlink);
1931 put_error_message("GM2CalcAmuTHDMMassBasis", "error",
1932 gm2calc_error_str(error));
1933 create_error_output();
1936 gm2calc_thdm_free(model);
1939/******************************************************************/
1941int main(int argc, char *argv[])
1943 gm2calc_sm_set_to_default(&sm);
1944 return MLMain(argc, argv);