BSMPT 3.1.3
BSMPT - Beyond the Standard Model Phase Transitions : A C++ package for the computation of the EWPT in BSM models
Loading...
Searching...
No Matches
ClassPotentialOrigin.h
Go to the documentation of this file.
1// Copyright (C) 2018 Philipp Basler and Margarete Mühlleitner
2// SPDX-FileCopyrightText: 2021 Philipp Basler, Margarete Mühlleitner and Jonas
3// Müller
4//
5// SPDX-License-Identifier: GPL-3.0-or-later
6
11#ifndef CLASSPOTENTIALORIGIN_H_
12#define CLASSPOTENTIALORIGIN_H_
13
14#include "Eigen/Eigenvalues"
18#include <BSMPT/utility/settings.h>
19#include <iostream>
20#include <vector>
21
22namespace BSMPT
23{
24
28const double C_PT = 0;
29
33const double C_threshold = 1e-4;
37static constexpr double C_CWcbFermion = 1.5;
41const double C_CWcbGB = 5.0 / 6.0;
45const double C_CWcbHiggs = 1.5;
46
53{
54public:
59
60protected:
64 bool UseTreeLevel = false;
65
69 double scale;
70
74 std::size_t nPar = 0;
78 std::size_t nParCT = 0;
79
83 std::vector<double> parStored;
84
88 std::vector<double> parCTStored;
89
93 ModelID::ModelIDs Model = ModelID::ModelIDs::NotSet;
97 std::size_t NNeutralHiggs = 0;
101 std::size_t NChargedHiggs = 0;
112 std::size_t NGauge = 4;
118 std::size_t NQuarks = 12;
119
126 std::size_t NColour = 3;
127
133 std::size_t NLepton = 9;
134
138 std::size_t nVEV = 0;
139
143 std::vector<double> VEVSymmetric;
147 std::vector<double> vevTree;
151 std::vector<double> vevTreeMin;
156 std::vector<double> TripleHiggsCorrectionsCW;
161 std::vector<std::vector<std::vector<double>>>
167 std::vector<std::vector<std::vector<double>>>
173 std::vector<std::vector<std::vector<double>>>
175
179 bool SetCurvatureDone = false;
184 bool CalcCouplingsdone = false;
190
195
199 bool UseIndexCol = false;
200
204 std::vector<double> HiggsVev;
208 std::vector<double> Curvature_Higgs_L1;
212 std::vector<std::vector<double>> Curvature_Higgs_L2;
216 std::vector<std::vector<std::vector<double>>> Curvature_Higgs_L3;
220 std::vector<std::vector<std::vector<std::vector<double>>>> Curvature_Higgs_L4;
224 std::vector<double> Curvature_Higgs_CT_L1;
228 std::vector<std::vector<double>> Curvature_Higgs_CT_L2;
232 std::vector<std::vector<std::vector<double>>> Curvature_Higgs_CT_L3;
236 std::vector<std::vector<std::vector<std::vector<double>>>>
241 std::vector<std::vector<std::vector<std::vector<double>>>>
246 std::vector<std::vector<std::vector<std::complex<double>>>>
251 std::vector<std::vector<std::complex<double>>> Curvature_Quark_F2;
255 std::vector<std::vector<std::vector<std::complex<double>>>>
260 std::vector<std::vector<std::complex<double>>> Curvature_Lepton_F2;
265 std::vector<double> MassSquaredHiggs;
270 std::vector<double> MassSquaredGauge;
275 std::vector<double> MassSquaredQuark;
280 std::vector<double> MassSquaredLepton;
285 std::vector<std::vector<double>> HiggsRotationMatrix;
290 std::vector<std::vector<std::vector<std::vector<double>>>>
296 std::vector<std::vector<std::vector<double>>> Couplings_Higgs_Triple;
301 std::vector<std::vector<std::vector<std::vector<double>>>>
307 std::vector<std::vector<std::vector<double>>> Couplings_Gauge_Higgs_21;
312 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
318 std::vector<std::vector<std::vector<std::complex<double>>>>
324 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
330 std::vector<std::vector<std::vector<std::complex<double>>>>
332
336 std::vector<std::vector<std::vector<double>>> LambdaGauge_3;
340 std::vector<std::vector<std::vector<double>>> LambdaHiggs_3;
345 std::vector<std::vector<std::vector<double>>> LambdaHiggs_3_CT;
350 std::vector<std::vector<std::vector<std::complex<double>>>> LambdaQuark_3;
355 std::vector<std::vector<std::vector<std::complex<double>>>> LambdaLepton_3;
361 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
368 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
370
375 std::vector<std::vector<double>> DebyeHiggs;
380 std::vector<std::vector<double>> DebyeGauge;
385 std::vector<std::size_t> VevOrder;
386
387public:
388 [[deprecated("Will call Class_Potential_Origin with GetSMConstants(). "
389 "Please use the "
390 "detailed overload "
391 "to ensure consistent SM constants through all "
392 "routines.")]] Class_Potential_Origin();
393
394 Class_Potential_Origin(const ISMConstants &smConstants);
395 virtual ~Class_Potential_Origin();
396
402 const auto &Get_Curvature_Higgs_L2() const { return Curvature_Higgs_L2; }
403
409 const auto &Get_Curvature_Higgs_L3() const { return Curvature_Higgs_L3; }
410
416 const auto &Get_Curvature_Higgs_L4() const { return Curvature_Higgs_L4; }
417
423 const auto &Get_Curvature_Gauge_G2H2() const { return Curvature_Gauge_G2H2; }
424
430 const auto &Get_Curvature_Lepton_F2H1() const
431 {
433 }
434
440 const auto &Get_Curvature_Lepton_F2() const { return Curvature_Lepton_F2; }
441
447 const auto &Get_Curvature_Quark_F2H1() const { return Curvature_Quark_F2H1; }
448
454 const auto &Get_Curvature_Quark_F2() const { return Curvature_Quark_F2; }
455
460 const std::vector<std::size_t> &Get_VevOrder() const { return VevOrder; }
461
466 double get_scale() const { return scale; }
467
472 void set_scale(double scale_new) { scale = scale_new; }
477 std::size_t get_nPar() const { return nPar; }
482 std::size_t get_nParCT() const { return nParCT; }
487 std::size_t get_nVEV() const { return nVEV; }
492 std::vector<double> get_vevTreeMin() const { return vevTreeMin; }
498 double get_vevTreeMin(const std::size_t &k) const { return vevTreeMin.at(k); }
503 std::vector<double> get_parStored() const { return parStored; }
508 std::vector<double> get_parCTStored() const { return parCTStored; }
513 std::size_t get_NGauge() const { return NGauge; }
518 std::size_t get_NQuarks() const { return NQuarks; }
523 std::size_t get_NHiggs() const { return NHiggs; }
528 std::size_t get_NChargedHiggs() const { return NChargedHiggs; }
533 std::size_t get_NNeutralHiggs() const { return NNeutralHiggs; }
538 std::size_t get_NLepton() const { return NLepton; }
543 ModelID::ModelIDs get_Model() const { return Model; }
544
549 const std::vector<std::vector<double>> &get_DebyeHiggs() const
550 {
551 return DebyeHiggs;
552 }
553
558 void set_InputLineNumber(int InputLineNumber_in)
559 {
560 InputLineNumber = InputLineNumber_in;
561 }
570 std::size_t j,
571 std::size_t k) const
572 {
573 return TripleHiggsCorrectionsTreePhysical.at(i).at(j).at(k);
574 }
583 std::size_t j,
584 std::size_t k) const
585 {
586 return TripleHiggsCorrectionsCTPhysical.at(i).at(j).at(k);
587 }
596 std::size_t j,
597 std::size_t k) const
598 {
599 return TripleHiggsCorrectionsCWPhysical.at(i).at(j).at(k);
600 }
601
605 void setUseIndexCol(std::string legend);
606
610 bool getUseIndexCol();
611
618 void sym2Dim(std::vector<std::vector<double>> &Tensor2Dim,
619 std::size_t Nk1,
620 std::size_t Nk2);
621
629 void sym3Dim(std::vector<std::vector<std::vector<double>>> &Tensor3Dim,
630 std::size_t Nk1,
631 std::size_t Nk2,
632 std::size_t Nk3);
633
642 void sym4Dim(
643 std::vector<std::vector<std::vector<std::vector<double>>>> &Tensor4Dim,
644 std::size_t Nk1,
645 std::size_t Nk2,
646 std::size_t Nk3,
647 std::size_t Nk4);
648
652 void initVectors();
653
663 virtual double VEff(const std::vector<double> &v,
664 double Temp = 0,
665 int diff = 0,
666 int Order = 1) const;
676 double VTree(const std::vector<double> &v,
677 int diff = 0,
678 bool ForceExplicitCalculation = false) const;
688 double CounterTerm(const std::vector<double> &v,
689 int diff = 0,
690 bool ForceExplicitCalculation = false) const;
701 double V1Loop(const std::vector<double> &v, double Temp, int diff) const;
702
707 virtual double EWSBVEV(const std::vector<double> &v) const;
708
715 virtual void SetEWVEVZero(std::vector<double> &sol) const;
716
720 virtual void ReadAndSet(const std::string &linestr,
721 std::vector<double> &par) = 0;
726 virtual std::vector<std::string> addLegendCT() const = 0;
732 virtual std::vector<std::string> addLegendTemp() const = 0;
737 virtual std::vector<std::string> addLegendTripleCouplings() const = 0;
742 virtual std::vector<std::string> addLegendVEV() const = 0;
743
749 virtual void set_gen(const std::vector<double> &par) = 0;
756 virtual void set_CT_Pot_Par(const std::vector<double> &par) = 0;
761 virtual void write() const = 0;
762
763 void set_All(const std::vector<double> &par,
764 const std::vector<double> &parCT);
765
771 virtual void SetCurvatureArrays() = 0;
780 std::vector<double> WeinbergFirstDerivative() const;
785 std::vector<double> WeinbergSecondDerivative() const;
790 Eigen::MatrixXd WeinbergSecondDerivativeAsMatrixXd() const;
795 std::vector<double> WeinbergThirdDerivative() const;
800 std::vector<double> WeinbergForthDerivative() const;
801
810 void CalculateDebye(bool forceCalculation = false);
811
816 void CalculateDebyeGauge();
821 void Prepare_Triple();
822
829 virtual bool CalculateDebyeSimplified() = 0;
830
838
843 virtual double VTreeSimplified(const std::vector<double> &v) const = 0;
849 bool UseVTreeSimplified = false;
854 virtual double VCounterSimplified(const std::vector<double> &v) const = 0;
861
873 std::vector<double> HiggsMassesSquared(const std::vector<double> &v,
874 const double &Temp = 0,
875 const int &diff = 0) const;
876
887 Eigen::MatrixXd HiggsMassMatrix(const std::vector<double> &v,
888 double Temp = 0,
889 int diff = 0) const;
890
902 std::vector<double> GaugeMassesSquared(const std::vector<double> &v,
903 const double &Temp = 0,
904 const int &diff = 0) const;
915 std::vector<double> QuarkMassesSquared(const std::vector<double> &v,
916 const int &diff = 0) const;
926 std::vector<double> LeptonMassesSquared(const std::vector<double> &v,
927 const int &diff = 0) const;
928
937 std::vector<std::complex<double>>
938 QuarkMasses(const std::vector<double> &v) const;
947 Eigen::MatrixXcd QuarkMassMatrix(const std::vector<double> &v) const;
956 std::vector<std::complex<double>>
957 LeptonMasses(const std::vector<double> &v) const;
966 Eigen::MatrixXcd LeptonMassMatrix(const std::vector<double> &v) const;
967
974 virtual void TripleHiggsCouplings() = 0;
979 double fbase(double MassSquaredA, double MassSquaredB) const;
984 double
985 fbaseTri(double MassSquaredA, double MassSquaredB, double MassSquaredC) const;
991 double fbaseFour(double MassSquaredA,
992 double MassSquaredB,
993 double MassSquaredC,
994 double MassSquaredD) const;
995
1000 virtual std::vector<double> calc_CT() const = 0;
1001
1007 double CWTerm(double MassSquared, double cb, int diff = 0) const;
1012 double FCW(double MassSquared) const;
1027 double boson(double MassSquared, double Temp, double cb, int diff = 0) const;
1032 [[deprecated(
1033 "Use this only if you want to calculate the effective potential with the "
1034 "exact same routine used in BSMPT v1.x. Use boson otherwise.")]] double
1035 boson_legacy(double MassSquared, double Temp, double cb, int diff = 0) const;
1046 double fermion(double MassSquared, double Temp, int diff = 0) const;
1051 [[deprecated(
1052 "Use this only if you want to calculate the effective potential with the "
1053 "exact same routine used in BSMPT v1.x. Use fermion otherwise.")]] double
1054 fermion_legacy(double Mass, double Temp, int diff = 0) const;
1059 [[deprecated("Use this only if you want to calculate the effective potential "
1060 "with the exact same routine used in BSMPT v1.x. "
1061 "Otherwise use ThermalFunctions::JInterpolatedHigh.")]] double
1062 Vl(double MassSquared, double Temp, int n, int diff = 0) const;
1071 [[deprecated(
1072 "Use this only if you want to calculate the effective potential with the "
1073 "exact same routine used in BSMPT v1.x. "
1074 "Otherwise use ThermalFunctions::JfermionInterpolatedLow.")]] double
1075 Vsf(double MassSquared, double Temp, int n, int diff = 0) const;
1076
1081 [[deprecated(
1082 "Use this only if you want to calculate the effective potential with the "
1083 "exact same routine used in BSMPT v1.x. "
1084 "Otherwise use ThermalFunctions::JbosonInterpolatedLow.")]] double
1085 Vsb(double MassSquaredIn, double Temp, int n, int diff = 0) const;
1086
1091 void resetbools();
1092
1102 std::vector<double>
1103 MinimizeOrderVEV(const std::vector<double> &VEVminimizer) const;
1104
1113 std::vector<double>
1114 FirstDerivativeOfEigenvalues(const Eigen::Ref<Eigen::MatrixXcd> M,
1115 const Eigen::Ref<Eigen::MatrixXcd> MDiff) const;
1129 std::vector<double> SecondDerivativeOfEigenvaluesNonRepeated(
1130 const Eigen::Ref<Eigen::MatrixXd> M,
1131 const Eigen::Ref<Eigen::MatrixXd> MDiffX,
1132 const Eigen::Ref<Eigen::MatrixXd> MDiffY,
1133 const Eigen::Ref<Eigen::MatrixXd> MDiffXY) const;
1134
1139 bool CheckNLOVEV(const std::vector<double> &v) const;
1140
1144 virtual void Debugging(const std::vector<double> &input,
1145 std::vector<double> &output) const = 0;
1146
1151 std::vector<std::vector<double>> SignSymmetries;
1156 void FindSignSymmetries();
1157
1161 void SetUseTreeLevel(bool val);
1162
1168 std::pair<std::vector<double>, std::vector<double>>
1169 initModel(std::string linestr);
1170
1179 std::vector<double> initModel(const std::vector<double> &par);
1180
1186 std::vector<double> resetScale(const double &newScale);
1187
1193 double CalculateRatioAlpha(const std::vector<double> &vev_symmetric,
1194 const std::vector<double> &vev_broken,
1195 const double &Temp) const;
1196
1201 Eigen::VectorXd NablaVCT(const std::vector<double> &v) const;
1206 Eigen::MatrixXd HessianCT(const std::vector<double> &v) const;
1207
1213 virtual std::vector<double> GetCTIdentities() const;
1214};
1215
1216} // namespace BSMPT
1217#endif /* CLASSPOTENTIALORIGIN_H_ */
The Class_Potential_Origin class Base class for all models. This class contains all numerical calcula...
Definition ClassPotentialOrigin.h:53
const auto & Get_Curvature_Quark_F2() const
Get_Curvature_Quark_F2 allows const ref access to the Y^{IJ} Tensor for Quarks.
Definition ClassPotentialOrigin.h:454
ModelID::ModelIDs get_Model() const
get_Model
Definition ClassPotentialOrigin.h:543
std::size_t nVEV
Definition ClassPotentialOrigin.h:138
const auto & Get_Curvature_Higgs_L2() const
Get_Curvature_Higgs_L2 allows const ref access to the L_{(S)}^{ij} Tensor.
Definition ClassPotentialOrigin.h:402
std::vector< double > WeinbergSecondDerivative() const
Definition ClassPotentialOrigin.cpp:1640
std::vector< double > MassSquaredGauge
MassSquaredGauge Stores the masses of the gauge Bosons calculated in CalculatePhysicalCouplings.
Definition ClassPotentialOrigin.h:270
std::vector< std::vector< double > > HiggsRotationMatrix
Definition ClassPotentialOrigin.h:285
std::vector< std::vector< double > > SignSymmetries
Definition ClassPotentialOrigin.h:1151
void CalculatePhysicalCouplings()
Definition ClassPotentialOrigin.cpp:962
void set_scale(double scale_new)
set_scale sets the MSBar renormalisation scale to scale_new
Definition ClassPotentialOrigin.h:472
bool CheckNLOVEV(const std::vector< double > &v) const
Definition ClassPotentialOrigin.cpp:3481
virtual double VTreeSimplified(const std::vector< double > &v) const =0
std::vector< double > Curvature_Higgs_L1
Definition ClassPotentialOrigin.h:208
double get_TripleHiggsCorrectionsCTPhysical(std::size_t i, std::size_t j, std::size_t k) const
get_TripleHiggsCorrectionsCTPhysical
Definition ClassPotentialOrigin.h:582
std::vector< double > get_parStored() const
get_parStored
Definition ClassPotentialOrigin.h:503
std::vector< double > vevTreeMin
Definition ClassPotentialOrigin.h:151
std::vector< double > get_parCTStored() const
get_parCTStored
Definition ClassPotentialOrigin.h:508
std::vector< std::vector< double > > DebyeGauge
DebyeGauge Stores the debye corrections to the mass matrix of the gauge bosons.
Definition ClassPotentialOrigin.h:380
std::vector< double > GaugeMassesSquared(const std::vector< double > &v, const double &Temp=0, const int &diff=0) const
Definition ClassPotentialOrigin.cpp:2510
double CWTerm(double MassSquared, double cb, int diff=0) const
Definition ClassPotentialOrigin.cpp:96
std::vector< std::vector< std::vector< std::vector< double > > > > Couplings_Higgs_Quartic
Couplings_Higgs_Quartic Stores the quartic Higgs couplings in the mass base.
Definition ClassPotentialOrigin.h:291
std::vector< double > FirstDerivativeOfEigenvalues(const Eigen::Ref< Eigen::MatrixXcd > M, const Eigen::Ref< Eigen::MatrixXcd > MDiff) const
Definition ClassPotentialOrigin.cpp:170
bool UseVTreeSimplified
UseVTreeSimplified Decides wether VTreeSimplified will be used or not. VTreeSimplified returns 0 if U...
Definition ClassPotentialOrigin.h:849
std::vector< std::vector< std::vector< double > > > TripleHiggsCorrectionsCTPhysical
Definition ClassPotentialOrigin.h:174
double VTree(const std::vector< double > &v, int diff=0, bool ForceExplicitCalculation=false) const
Definition ClassPotentialOrigin.cpp:2826
virtual void ReadAndSet(const std::string &linestr, std::vector< double > &par)=0
std::vector< std::vector< std::vector< double > > > Curvature_Higgs_L3
Definition ClassPotentialOrigin.h:216
std::vector< std::vector< std::vector< std::vector< std::complex< double > > > > > Couplings_Lepton_Higgs_22
Couplings_Lepton_Higgs_22 Stores the couplings between two leptons and two Higgs bosons in the mass b...
Definition ClassPotentialOrigin.h:325
std::vector< double > LeptonMassesSquared(const std::vector< double > &v, const int &diff=0) const
Definition ClassPotentialOrigin.cpp:2729
std::vector< double > WeinbergThirdDerivative() const
Definition ClassPotentialOrigin.cpp:1656
std::size_t get_nVEV() const
get_nVEV
Definition ClassPotentialOrigin.h:487
void resetbools()
Definition ClassPotentialOrigin.cpp:3472
std::vector< std::vector< std::vector< std::complex< double > > > > Couplings_Quark_Higgs_21
Couplings_Quark_Higgs_21 Stores the couplings between two quarks and one Higgs boson in the mass base...
Definition ClassPotentialOrigin.h:319
std::vector< double > SecondDerivativeOfEigenvaluesNonRepeated(const Eigen::Ref< Eigen::MatrixXd > M, const Eigen::Ref< Eigen::MatrixXd > MDiffX, const Eigen::Ref< Eigen::MatrixXd > MDiffY, const Eigen::Ref< Eigen::MatrixXd > MDiffXY) const
Definition ClassPotentialOrigin.cpp:887
std::size_t nPar
Definition ClassPotentialOrigin.h:74
virtual bool CalculateDebyeSimplified()=0
void sym2Dim(std::vector< std::vector< double > > &Tensor2Dim, std::size_t Nk1, std::size_t Nk2)
sym2Dim Symmetrize scalar 2-dim tensor
Definition ClassPotentialOrigin.cpp:3392
std::vector< std::vector< std::vector< double > > > Couplings_Gauge_Higgs_21
Couplings_Gauge_Higgs_21 Stores the coupling between two gauge and one Higgs boson in the mass base.
Definition ClassPotentialOrigin.h:307
std::vector< double > WeinbergForthDerivative() const
Definition ClassPotentialOrigin.cpp:1902
std::vector< std::complex< double > > QuarkMasses(const std::vector< double > &v) const
Definition ClassPotentialOrigin.cpp:3687
const ISMConstants SMConstants
SMConstants The SM constants used by the model.
Definition ClassPotentialOrigin.h:58
std::vector< std::vector< double > > Curvature_Higgs_CT_L2
Definition ClassPotentialOrigin.h:228
std::vector< double > Curvature_Higgs_CT_L1
Definition ClassPotentialOrigin.h:224
Eigen::MatrixXd HiggsMassMatrix(const std::vector< double > &v, double Temp=0, int diff=0) const
HiggsMassMatrix calculates the Higgs mass matrix.
Definition ClassPotentialOrigin.cpp:2324
std::vector< double > vevTree
Definition ClassPotentialOrigin.h:147
virtual std::vector< double > calc_CT() const =0
std::vector< std::vector< std::vector< std::vector< double > > > > Curvature_Gauge_G2H2
Definition ClassPotentialOrigin.h:242
std::vector< std::vector< std::vector< std::vector< std::complex< double > > > > > LambdaQuark_4
LambdaQuark_4 Stores the Lambda_{(F)}^{IJkm} tensor for quarks , describing the derivative of Lambda_...
Definition ClassPotentialOrigin.h:362
std::vector< std::vector< std::vector< double > > > TripleHiggsCorrectionsCWPhysical
Definition ClassPotentialOrigin.h:162
void CalculateDebyeGauge()
Definition ClassPotentialOrigin.cpp:3292
double get_TripleHiggsCorrectionsTreePhysical(std::size_t i, std::size_t j, std::size_t k) const
get_TripleHiggsCorrectionsTreePhysical
Definition ClassPotentialOrigin.h:569
std::vector< std::vector< std::vector< std::vector< std::complex< double > > > > > Couplings_Quark_Higgs_22
Couplings_Quark_Higgs_22 Stores the couplings between two quarks and two Higgs bosons in the mass bas...
Definition ClassPotentialOrigin.h:313
void CalculateDebye(bool forceCalculation=false)
Definition ClassPotentialOrigin.cpp:3223
std::vector< std::vector< std::vector< std::vector< double > > > > Curvature_Higgs_L4
Definition ClassPotentialOrigin.h:220
std::vector< std::vector< std::vector< double > > > Curvature_Higgs_CT_L3
Definition ClassPotentialOrigin.h:232
Eigen::MatrixXd WeinbergSecondDerivativeAsMatrixXd() const
Definition ClassPotentialOrigin.cpp:1498
double fbaseFour(double MassSquaredA, double MassSquaredB, double MassSquaredC, double MassSquaredD) const
Definition ClassPotentialOrigin.cpp:389
virtual double VCounterSimplified(const std::vector< double > &v) const =0
std::vector< std::size_t > VevOrder
VevOrder Stores the matching order used in MinimizeOrderVEV, set in the constructor of the model.
Definition ClassPotentialOrigin.h:385
std::size_t nParCT
Definition ClassPotentialOrigin.h:78
void setUseIndexCol(std::string legend)
Definition ClassPotentialOrigin.cpp:3534
bool CalculatedTripleCopulings
CalculatedTripleCopulings Used to check if TripleHiggsCouplings has already been called.
Definition ClassPotentialOrigin.h:189
double fbase(double MassSquaredA, double MassSquaredB) const
Definition ClassPotentialOrigin.cpp:859
std::vector< std::vector< std::vector< double > > > LambdaHiggs_3_CT
LambdaHiggs_3 Stores the Lambda_{(S)}^{ijk} tensor for the counterterm parameters.
Definition ClassPotentialOrigin.h:345
std::vector< std::vector< std::vector< std::complex< double > > > > Curvature_Lepton_F2H1
Definition ClassPotentialOrigin.h:256
Eigen::MatrixXcd QuarkMassMatrix(const std::vector< double > &v) const
QuarkMassMatrix calculates the Mass Matrix for the Quarks of the form $ M^{IJ} = Y^{IJ} + Y^{IJk} v_k...
Definition ClassPotentialOrigin.cpp:3633
std::vector< std::vector< std::vector< double > > > TripleHiggsCorrectionsTreePhysical
Definition ClassPotentialOrigin.h:168
double boson_legacy(double MassSquared, double Temp, double cb, int diff=0) const
Definition ClassPotentialOrigin_deprecated.cpp:221
void sym4Dim(std::vector< std::vector< std::vector< std::vector< double > > > > &Tensor4Dim, std::size_t Nk1, std::size_t Nk2, std::size_t Nk3, std::size_t Nk4)
sym4Dim Symmetrize scalar 4-dim tensor
Definition ClassPotentialOrigin.cpp:3428
const std::vector< std::size_t > & Get_VevOrder() const
Get_VevOrder allows const ref access to the VevOrder vector.
Definition ClassPotentialOrigin.h:460
virtual void set_CT_Pot_Par(const std::vector< double > &par)=0
virtual void TripleHiggsCouplings()=0
double Vl(double MassSquared, double Temp, int n, int diff=0) const
Definition ClassPotentialOrigin_deprecated.cpp:31
virtual void Debugging(const std::vector< double > &input, std::vector< double > &output) const =0
std::size_t NChargedHiggs
Definition ClassPotentialOrigin.h:101
bool getUseIndexCol()
Definition ClassPotentialOrigin.cpp:3539
double get_TripleHiggsCorrectionsCWPhysical(std::size_t i, std::size_t j, std::size_t k) const
get_TripleHiggsCorrectionsCWPhysical
Definition ClassPotentialOrigin.h:595
virtual void write() const =0
std::vector< std::vector< std::vector< std::complex< double > > > > Curvature_Quark_F2H1
Definition ClassPotentialOrigin.h:247
std::vector< std::vector< std::vector< std::complex< double > > > > LambdaLepton_3
LambdaLepton_3 Stores the Lambda_{(F)}^{IJk} tensor for leptons , describing the derivative of Lambda...
Definition ClassPotentialOrigin.h:355
double scale
Definition ClassPotentialOrigin.h:69
std::vector< double > parStored
Definition ClassPotentialOrigin.h:83
bool CalcCouplingsdone
CalcCouplingsdone Used to check if CalculatePhysicalCouplings has already been called.
Definition ClassPotentialOrigin.h:184
std::size_t NColour
NColour Number of colours of the quarks. Do not change this in the current version as we only investi...
Definition ClassPotentialOrigin.h:126
void Prepare_Triple()
Definition ClassPotentialOrigin.cpp:61
double Vsb(double MassSquaredIn, double Temp, int n, int diff=0) const
Definition ClassPotentialOrigin_deprecated.cpp:149
std::vector< double > MinimizeOrderVEV(const std::vector< double > &VEVminimizer) const
Definition ClassPotentialOrigin.cpp:3800
std::size_t get_NChargedHiggs() const
get_NChargedHiggs
Definition ClassPotentialOrigin.h:528
std::vector< std::vector< std::complex< double > > > Curvature_Lepton_F2
Definition ClassPotentialOrigin.h:260
ModelID::ModelIDs Model
Definition ClassPotentialOrigin.h:93
virtual std::vector< std::string > addLegendTemp() const =0
std::vector< std::vector< std::vector< double > > > LambdaGauge_3
LambdaGauge_3 Stores the Lambda_{(G)}^{abi} tensor.
Definition ClassPotentialOrigin.h:336
const std::vector< std::vector< double > > & get_DebyeHiggs() const
get_DebyeHiggs get the Debye corrections to the Higgs mass matrix
Definition ClassPotentialOrigin.h:549
std::size_t NNeutralHiggs
Definition ClassPotentialOrigin.h:97
void initVectors()
Definition ClassPotentialOrigin.cpp:3332
std::size_t get_NNeutralHiggs() const
get_NNeutralHiggs
Definition ClassPotentialOrigin.h:533
void sym3Dim(std::vector< std::vector< std::vector< double > > > &Tensor3Dim, std::size_t Nk1, std::size_t Nk2, std::size_t Nk3)
sym3Dim Symmetrize scalar 3-dim tensor
Definition ClassPotentialOrigin.cpp:3406
bool SetCurvatureDone
SetCurvatureDone Used to check if the tensors are set.
Definition ClassPotentialOrigin.h:179
double FCW(double MassSquared) const
Definition ClassPotentialOrigin.cpp:80
int InputLineNumber
InputLineNumber Used for the error message in fbaseTri.
Definition ClassPotentialOrigin.h:194
Eigen::VectorXd NablaVCT(const std::vector< double > &v) const
NablaVCT.
Definition ClassPotentialOrigin.cpp:3820
std::vector< std::vector< std::vector< double > > > Couplings_Higgs_Triple
Couplings_Higgs_Triple Stores the triple Higgs couplings in the mass base.
Definition ClassPotentialOrigin.h:296
virtual std::vector< std::string > addLegendCT() const =0
std::vector< double > MassSquaredLepton
MassSquaredLepton Stores the masses of the leptons calculated in CalculatePhysicalCouplings.
Definition ClassPotentialOrigin.h:280
std::vector< double > MassSquaredHiggs
MassSquaredHiggs Stores the masses of the Higgs Bosons calculated in CalculatePhysicalCouplings.
Definition ClassPotentialOrigin.h:265
double get_vevTreeMin(const std::size_t &k) const
get_vevTreeMin
Definition ClassPotentialOrigin.h:498
virtual double EWSBVEV(const std::vector< double > &v) const
Definition ClassPotentialOrigin.cpp:3495
const auto & Get_Curvature_Quark_F2H1() const
Get_Curvature_Quark_F2H1 allows const ref access to the Y^{IJk} Tensor for Quarks.
Definition ClassPotentialOrigin.h:447
bool UseVCounterSimplified
UseVCounterSimplified Decides wether VCounterSimplified will be used or not. VCounterSimplified retur...
Definition ClassPotentialOrigin.h:860
const auto & Get_Curvature_Lepton_F2H1() const
Get_Curvature_Lepton_F2H1 allows const ref access to the Y^{IJk} Tensor for Leptons.
Definition ClassPotentialOrigin.h:430
const auto & Get_Curvature_Lepton_F2() const
Get_Curvature_Lepton_F2 allows const ref access to the Y^{IJ} Tensor for Leptons.
Definition ClassPotentialOrigin.h:440
bool UseTreeLevel
UseTreeLevel Enforces VEff to only use the tree-level potential.
Definition ClassPotentialOrigin.h:64
double fermion(double MassSquared, double Temp, int diff=0) const
Calculation of the fermionic std::size_tegral + Coleman-Weinberg potential without taking d....
Definition ClassPotentialOrigin.cpp:143
std::size_t NQuarks
Definition ClassPotentialOrigin.h:118
double get_scale() const
get_scale
Definition ClassPotentialOrigin.h:466
std::vector< double > WeinbergFirstDerivative() const
Definition ClassPotentialOrigin.cpp:1403
std::vector< std::vector< std::vector< std::vector< double > > > > Couplings_Gauge_Higgs_22
Couplings_Gauge_Higgs_22 Stores the couplings between two Higgs and two gauge bosons in the mass base...
Definition ClassPotentialOrigin.h:302
void SetUseTreeLevel(bool val)
Definition ClassPotentialOrigin.cpp:3580
virtual std::vector< std::string > addLegendTripleCouplings() const =0
std::vector< std::vector< double > > DebyeHiggs
DebyeHiggs Stores the debye corrections to the mass matrix of the Higgs bosons.
Definition ClassPotentialOrigin.h:375
std::vector< double > parCTStored
Definition ClassPotentialOrigin.h:88
std::size_t get_nPar() const
get_nPar
Definition ClassPotentialOrigin.h:477
double boson(double MassSquared, double Temp, double cb, int diff=0) const
Calculation of the bosonic std::size_tegral + Coleman-Weinberg potential without taking d....
Definition ClassPotentialOrigin.cpp:112
virtual std::vector< std::string > addLegendVEV() const =0
std::vector< std::vector< std::vector< double > > > LambdaHiggs_3
LambdaHiggs_3 Stores the Lambda_{(S)}^{ijk} tensor.
Definition ClassPotentialOrigin.h:340
const auto & Get_Curvature_Gauge_G2H2() const
Get_Curvature_Gauge_G2H2 allows const ref access to the G^{abij} Tensor.
Definition ClassPotentialOrigin.h:423
std::size_t NLepton
Definition ClassPotentialOrigin.h:133
virtual void set_gen(const std::vector< double > &par)=0
std::size_t NHiggs
Definition ClassPotentialOrigin.h:106
virtual void SetEWVEVZero(std::vector< double > &sol) const
SetEWVEVZero Set all VEV directions in sol-vector to zero that contibute to EW VEV.
Definition ClassPotentialOrigin.cpp:3516
virtual void SetCurvatureArrays()=0
std::vector< double > QuarkMassesSquared(const std::vector< double > &v, const int &diff=0) const
Definition ClassPotentialOrigin.cpp:2631
std::vector< std::complex< double > > LeptonMasses(const std::vector< double > &v) const
Definition ClassPotentialOrigin.cpp:3761
double V1Loop(const std::vector< double > &v, double Temp, int diff) const
Definition ClassPotentialOrigin.cpp:2991
std::vector< std::vector< std::vector< std::vector< double > > > > Curvature_Higgs_CT_L4
Definition ClassPotentialOrigin.h:237
std::vector< std::vector< std::vector< std::complex< double > > > > LambdaQuark_3
LambdaQuark_3 Stores the Lambda_{(F)}^{IJk} tensor for quarks , describing the derivative of Lambda_{...
Definition ClassPotentialOrigin.h:350
std::vector< std::vector< std::vector< std::complex< double > > > > Couplings_Lepton_Higgs_21
Couplings_Quark_Higgs_21 Stores the couplings between two leptons and one Higgs boson in the mass bas...
Definition ClassPotentialOrigin.h:331
std::vector< double > HiggsVev
Definition ClassPotentialOrigin.h:204
std::size_t get_nParCT() const
get_nParCT
Definition ClassPotentialOrigin.h:482
virtual std::vector< double > GetCTIdentities() const
GetCTIdentities.
Definition ClassPotentialOrigin.cpp:3866
std::vector< double > HiggsMassesSquared(const std::vector< double > &v, const double &Temp=0, const int &diff=0) const
Definition ClassPotentialOrigin.cpp:2421
void set_All(const std::vector< double > &par, const std::vector< double > &parCT)
Definition ClassPotentialOrigin.cpp:50
std::vector< double > MassSquaredQuark
MassSquaredQuark Stores the masses of the quarks calculated in CalculatePhysicalCouplings.
Definition ClassPotentialOrigin.h:275
virtual double VEff(const std::vector< double > &v, double Temp=0, int diff=0, int Order=1) const
Definition ClassPotentialOrigin.cpp:2948
double CalculateRatioAlpha(const std::vector< double > &vev_symmetric, const std::vector< double > &vev_broken, const double &Temp) const
Definition ClassPotentialOrigin.cpp:3780
double fbaseTri(double MassSquaredA, double MassSquaredB, double MassSquaredC) const
Definition ClassPotentialOrigin.cpp:288
std::vector< double > get_vevTreeMin() const
get_vevTreeMin
Definition ClassPotentialOrigin.h:492
std::size_t get_NQuarks() const
get_NQuarks
Definition ClassPotentialOrigin.h:518
const auto & Get_Curvature_Higgs_L3() const
Get_Curvature_Higgs_L3 allows const ref access to the L_{(S)}^{ijk} Tensor.
Definition ClassPotentialOrigin.h:409
virtual bool CalculateDebyeGaugeSimplified()=0
std::size_t get_NHiggs() const
get_NHiggs
Definition ClassPotentialOrigin.h:523
std::vector< double > VEVSymmetric
Definition ClassPotentialOrigin.h:143
double Vsf(double MassSquared, double Temp, int n, int diff=0) const
Definition ClassPotentialOrigin_deprecated.cpp:83
std::pair< std::vector< double >, std::vector< double > > initModel(std::string linestr)
Definition ClassPotentialOrigin.cpp:3586
Eigen::MatrixXd HessianCT(const std::vector< double > &v) const
HessianWeinberg.
Definition ClassPotentialOrigin.cpp:3844
std::vector< std::vector< double > > Curvature_Higgs_L2
Definition ClassPotentialOrigin.h:212
std::vector< std::vector< std::vector< std::vector< std::complex< double > > > > > LambdaLepton_4
LambdaLepton_4 Stores the Lambda_{(F)}^{IJkm} tensor for leptons , describing the derivative of Lambd...
Definition ClassPotentialOrigin.h:369
std::size_t NGauge
Definition ClassPotentialOrigin.h:112
const auto & Get_Curvature_Higgs_L4() const
Get_Curvature_Higgs_L4 allows const ref access to the L_{(S)}^{ijkl} Tensor.
Definition ClassPotentialOrigin.h:416
std::size_t get_NLepton() const
get_NLepton
Definition ClassPotentialOrigin.h:538
std::size_t get_NGauge() const
get_NGauge
Definition ClassPotentialOrigin.h:513
std::vector< std::vector< std::complex< double > > > Curvature_Quark_F2
Definition ClassPotentialOrigin.h:251
bool UseIndexCol
Definition ClassPotentialOrigin.h:199
std::vector< double > resetScale(const double &newScale)
resetScale changes the MSBar scale to newScale
Definition ClassPotentialOrigin.cpp:3621
Eigen::MatrixXcd LeptonMassMatrix(const std::vector< double > &v) const
LeptonMassMatrix calculates the Mass Matrix for the Leptons of the form $ M^{IJ} = Y^{IJ} + Y^{IJk} v...
Definition ClassPotentialOrigin.cpp:3709
double CounterTerm(const std::vector< double > &v, int diff=0, bool ForceExplicitCalculation=false) const
Definition ClassPotentialOrigin.cpp:2892
double fermion_legacy(double Mass, double Temp, int diff=0) const
Definition ClassPotentialOrigin_deprecated.cpp:316
std::vector< double > TripleHiggsCorrectionsCW
Definition ClassPotentialOrigin.h:156
void FindSignSymmetries()
FindSignSymmetries checks for all possible sign changes in the VEV components and checks for all poss...
Definition ClassPotentialOrigin.cpp:3544
void set_InputLineNumber(int InputLineNumber_in)
set_InputLineNumber
Definition ClassPotentialOrigin.h:558
This classes calculates the Bounce action of the potential with a set temperature.
Definition CalculateEtaInterface.h:24
const double C_PT
C_PT Lower threshold to stop the EWPT calculation.
Definition ClassPotentialOrigin.h:28
const double C_CWcbHiggs
C_CWcbHiggs constant used in the CW potential for Higgs bosons.
Definition ClassPotentialOrigin.h:45
const double C_CWcbGB
C_CWcbGB constant used in the CW potential for gauge bosons.
Definition ClassPotentialOrigin.h:41
const double C_threshold
C_threshold threshold to check if a mass is numerically zero.
Definition ClassPotentialOrigin.h:33
The ISMConstants struct containing all necessary SM constants.
Definition SMparam.h:24
Definition transition_tracer.h:157