BSMPT 3.0.7
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 <iostream>
19#include <vector>
20
21namespace BSMPT
22{
23
27const bool C_UseParwani = false;
28
32const double C_PT = 0;
33
37const double C_threshold = 1e-4;
41static constexpr double C_CWcbFermion = 1.5;
45const double C_CWcbGB = 5.0 / 6.0;
49const double C_CWcbHiggs = 1.5;
50
57{
58public:
63
64protected:
68 bool UseTreeLevel = false;
69
73 double scale;
74
78 std::size_t nPar = 0;
82 std::size_t nParCT = 0;
83
87 std::vector<double> parStored;
88
92 std::vector<double> parCTStored;
93
97 ModelID::ModelIDs Model = ModelID::ModelIDs::NotSet;
101 std::size_t NNeutralHiggs = 0;
105 std::size_t NChargedHiggs = 0;
116 std::size_t NGauge = 4;
122 std::size_t NQuarks = 12;
123
130 std::size_t NColour = 3;
131
137 std::size_t NLepton = 9;
138
142 std::size_t nVEV = 0;
143
147 std::vector<double> VEVSymmetric;
151 std::vector<double> vevTree;
155 std::vector<double> vevTreeMin;
160 std::vector<double> TripleHiggsCorrectionsCW;
165 std::vector<std::vector<std::vector<double>>>
171 std::vector<std::vector<std::vector<double>>>
177 std::vector<std::vector<std::vector<double>>>
179
183 bool SetCurvatureDone = false;
188 bool CalcCouplingsdone = false;
194
199
203 bool UseIndexCol = false;
204
208 std::vector<double> HiggsVev;
212 std::vector<double> Curvature_Higgs_L1;
216 std::vector<std::vector<double>> Curvature_Higgs_L2;
220 std::vector<std::vector<std::vector<double>>> Curvature_Higgs_L3;
224 std::vector<std::vector<std::vector<std::vector<double>>>> Curvature_Higgs_L4;
228 std::vector<double> Curvature_Higgs_CT_L1;
232 std::vector<std::vector<double>> Curvature_Higgs_CT_L2;
236 std::vector<std::vector<std::vector<double>>> Curvature_Higgs_CT_L3;
240 std::vector<std::vector<std::vector<std::vector<double>>>>
245 std::vector<std::vector<std::vector<std::vector<double>>>>
250 std::vector<std::vector<std::vector<std::complex<double>>>>
255 std::vector<std::vector<std::complex<double>>> Curvature_Quark_F2;
259 std::vector<std::vector<std::vector<std::complex<double>>>>
264 std::vector<std::vector<std::complex<double>>> Curvature_Lepton_F2;
269 std::vector<double> MassSquaredHiggs;
274 std::vector<double> MassSquaredGauge;
279 std::vector<double> MassSquaredQuark;
284 std::vector<double> MassSquaredLepton;
289 std::vector<std::vector<double>> HiggsRotationMatrix;
294 std::vector<std::vector<std::vector<std::vector<double>>>>
300 std::vector<std::vector<std::vector<double>>> Couplings_Higgs_Triple;
305 std::vector<std::vector<std::vector<std::vector<double>>>>
311 std::vector<std::vector<std::vector<double>>> Couplings_Gauge_Higgs_21;
316 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
322 std::vector<std::vector<std::vector<std::complex<double>>>>
328 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
334 std::vector<std::vector<std::vector<std::complex<double>>>>
336
340 std::vector<std::vector<std::vector<double>>> LambdaGauge_3;
344 std::vector<std::vector<std::vector<double>>> LambdaHiggs_3;
349 std::vector<std::vector<std::vector<double>>> LambdaHiggs_3_CT;
354 std::vector<std::vector<std::vector<std::complex<double>>>> LambdaQuark_3;
359 std::vector<std::vector<std::vector<std::complex<double>>>> LambdaLepton_3;
365 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
372 std::vector<std::vector<std::vector<std::vector<std::complex<double>>>>>
374
379 std::vector<std::vector<double>> DebyeHiggs;
384 std::vector<std::vector<double>> DebyeGauge;
389 std::vector<std::size_t> VevOrder;
390
391public:
392 [[deprecated("Will call Class_Potential_Origin with GetSMConstants(). "
393 "Please use the "
394 "detailed overload "
395 "to ensure consistent SM constants through all "
396 "routines.")]] Class_Potential_Origin();
397
398 Class_Potential_Origin(const ISMConstants &smConstants);
399 virtual ~Class_Potential_Origin();
400
406 const auto &Get_Curvature_Higgs_L2() const { return Curvature_Higgs_L2; }
407
413 const auto &Get_Curvature_Higgs_L3() const { return Curvature_Higgs_L3; }
414
420 const auto &Get_Curvature_Higgs_L4() const { return Curvature_Higgs_L4; }
421
427 const auto &Get_Curvature_Gauge_G2H2() const { return Curvature_Gauge_G2H2; }
428
434 const auto &Get_Curvature_Lepton_F2H1() const
435 {
437 }
438
444 const auto &Get_Curvature_Lepton_F2() const { return Curvature_Lepton_F2; }
445
451 const auto &Get_Curvature_Quark_F2H1() const { return Curvature_Quark_F2H1; }
452
458 const auto &Get_Curvature_Quark_F2() const { return Curvature_Quark_F2; }
459
464 const std::vector<std::size_t> &Get_VevOrder() const { return VevOrder; }
465
470 double get_scale() const { return scale; }
471
476 void set_scale(double scale_new) { scale = scale_new; }
481 std::size_t get_nPar() const { return nPar; }
486 std::size_t get_nParCT() const { return nParCT; }
491 std::size_t get_nVEV() const { return nVEV; }
496 std::vector<double> get_vevTreeMin() const { return vevTreeMin; }
502 double get_vevTreeMin(const std::size_t &k) const { return vevTreeMin.at(k); }
507 std::vector<double> get_parStored() const { return parStored; }
512 std::vector<double> get_parCTStored() const { return parCTStored; }
517 std::size_t get_NGauge() const { return NGauge; }
522 std::size_t get_NQuarks() const { return NQuarks; }
527 std::size_t get_NHiggs() const { return NHiggs; }
532 std::size_t get_NChargedHiggs() const { return NChargedHiggs; }
537 std::size_t get_NNeutralHiggs() const { return NNeutralHiggs; }
542 std::size_t get_NLepton() const { return NLepton; }
548
553 const std::vector<std::vector<double>> &get_DebyeHiggs() const
554 {
555 return DebyeHiggs;
556 }
557
562 void set_InputLineNumber(int InputLineNumber_in)
563 {
564 InputLineNumber = InputLineNumber_in;
565 }
574 std::size_t j,
575 std::size_t k) const
576 {
577 return TripleHiggsCorrectionsTreePhysical.at(i).at(j).at(k);
578 }
587 std::size_t j,
588 std::size_t k) const
589 {
590 return TripleHiggsCorrectionsCTPhysical.at(i).at(j).at(k);
591 }
600 std::size_t j,
601 std::size_t k) const
602 {
603 return TripleHiggsCorrectionsCWPhysical.at(i).at(j).at(k);
604 }
605
609 void setUseIndexCol(std::string legend);
610
614 bool getUseIndexCol();
615
622 void sym2Dim(std::vector<std::vector<double>> &Tensor2Dim,
623 std::size_t Nk1,
624 std::size_t Nk2);
625
633 void sym3Dim(std::vector<std::vector<std::vector<double>>> &Tensor3Dim,
634 std::size_t Nk1,
635 std::size_t Nk2,
636 std::size_t Nk3);
637
646 void sym4Dim(
647 std::vector<std::vector<std::vector<std::vector<double>>>> &Tensor4Dim,
648 std::size_t Nk1,
649 std::size_t Nk2,
650 std::size_t Nk3,
651 std::size_t Nk4);
652
656 void initVectors();
657
667 virtual double VEff(const std::vector<double> &v,
668 double Temp = 0,
669 int diff = 0,
670 int Order = 1) const;
680 double VTree(const std::vector<double> &v,
681 int diff = 0,
682 bool ForceExplicitCalculation = false) const;
692 double CounterTerm(const std::vector<double> &v,
693 int diff = 0,
694 bool ForceExplicitCalculation = false) const;
705 double V1Loop(const std::vector<double> &v, double Temp, int diff) const;
706
711 virtual double EWSBVEV(const std::vector<double> &v) const;
712
719 virtual void SetEWVEVZero(std::vector<double> &sol) const;
720
724 virtual void ReadAndSet(const std::string &linestr,
725 std::vector<double> &par) = 0;
730 virtual std::vector<std::string> addLegendCT() const = 0;
736 virtual std::vector<std::string> addLegendTemp() const = 0;
741 virtual std::vector<std::string> addLegendTripleCouplings() const = 0;
746 virtual std::vector<std::string> addLegendVEV() const = 0;
747
753 virtual void set_gen(const std::vector<double> &par) = 0;
760 virtual void set_CT_Pot_Par(const std::vector<double> &par) = 0;
765 virtual void write() const = 0;
766
767 void set_All(const std::vector<double> &par,
768 const std::vector<double> &parCT);
769
775 virtual void SetCurvatureArrays() = 0;
784 std::vector<double> WeinbergFirstDerivative() const;
789 std::vector<double> WeinbergSecondDerivative() const;
794 Eigen::MatrixXd WeinbergSecondDerivativeAsMatrixXd() const;
799 std::vector<double> WeinbergThirdDerivative() const;
804 std::vector<double> WeinbergForthDerivative() const;
805
814 void CalculateDebye(bool forceCalculation = false);
815
820 void CalculateDebyeGauge();
825 void Prepare_Triple();
826
833 virtual bool CalculateDebyeSimplified() = 0;
834
842
847 virtual double VTreeSimplified(const std::vector<double> &v) const = 0;
853 bool UseVTreeSimplified = false;
858 virtual double VCounterSimplified(const std::vector<double> &v) const = 0;
865
877 std::vector<double> HiggsMassesSquared(const std::vector<double> &v,
878 const double &Temp = 0,
879 const int &diff = 0) const;
880
891 Eigen::MatrixXd HiggsMassMatrix(const std::vector<double> &v,
892 double Temp = 0,
893 int diff = 0) const;
894
906 std::vector<double> GaugeMassesSquared(const std::vector<double> &v,
907 const double &Temp = 0,
908 const int &diff = 0) const;
919 std::vector<double> QuarkMassesSquared(const std::vector<double> &v,
920 const int &diff = 0) const;
930 std::vector<double> LeptonMassesSquared(const std::vector<double> &v,
931 const int &diff = 0) const;
932
941 std::vector<std::complex<double>>
942 QuarkMasses(const std::vector<double> &v) const;
951 Eigen::MatrixXcd QuarkMassMatrix(const std::vector<double> &v) const;
960 std::vector<std::complex<double>>
961 LeptonMasses(const std::vector<double> &v) const;
970 Eigen::MatrixXcd LeptonMassMatrix(const std::vector<double> &v) const;
971
978 virtual void TripleHiggsCouplings() = 0;
983 double fbase(double MassSquaredA, double MassSquaredB) const;
988 double
989 fbaseTri(double MassSquaredA, double MassSquaredB, double MassSquaredC) const;
995 double fbaseFour(double MassSquaredA,
996 double MassSquaredB,
997 double MassSquaredC,
998 double MassSquaredD) const;
999
1004 virtual std::vector<double> calc_CT() const = 0;
1005
1011 double CWTerm(double MassSquared, double cb, int diff = 0) const;
1016 double FCW(double MassSquared) const;
1031 double boson(double MassSquared, double Temp, double cb, int diff = 0) const;
1036 [[deprecated(
1037 "Use this only if you want to calculate the effective potential with the "
1038 "exact same routine used in BSMPT v1.x. Use boson otherwise.")]] double
1039 boson_legacy(double MassSquared, double Temp, double cb, int diff = 0) const;
1050 double fermion(double MassSquared, double Temp, int diff = 0) const;
1055 [[deprecated(
1056 "Use this only if you want to calculate the effective potential with the "
1057 "exact same routine used in BSMPT v1.x. Use fermion otherwise.")]] double
1058 fermion_legacy(double Mass, double Temp, int diff = 0) const;
1063 [[deprecated("Use this only if you want to calculate the effective potential "
1064 "with the exact same routine used in BSMPT v1.x. "
1065 "Otherwise use ThermalFunctions::JInterpolatedHigh.")]] double
1066 Vl(double MassSquared, double Temp, int n, int diff = 0) const;
1075 [[deprecated(
1076 "Use this only if you want to calculate the effective potential with the "
1077 "exact same routine used in BSMPT v1.x. "
1078 "Otherwise use ThermalFunctions::JfermionInterpolatedLow.")]] double
1079 Vsf(double MassSquared, double Temp, int n, int diff = 0) const;
1080
1085 [[deprecated(
1086 "Use this only if you want to calculate the effective potential with the "
1087 "exact same routine used in BSMPT v1.x. "
1088 "Otherwise use ThermalFunctions::JbosonInterpolatedLow.")]] double
1089 Vsb(double MassSquaredIn, double Temp, int n, int diff = 0) const;
1090
1095 void resetbools();
1096
1106 std::vector<double>
1107 MinimizeOrderVEV(const std::vector<double> &VEVminimizer) const;
1108
1117 std::vector<double>
1118 FirstDerivativeOfEigenvalues(const Eigen::Ref<Eigen::MatrixXcd> M,
1119 const Eigen::Ref<Eigen::MatrixXcd> MDiff) const;
1133 std::vector<double> SecondDerivativeOfEigenvaluesNonRepeated(
1134 const Eigen::Ref<Eigen::MatrixXd> M,
1135 const Eigen::Ref<Eigen::MatrixXd> MDiffX,
1136 const Eigen::Ref<Eigen::MatrixXd> MDiffY,
1137 const Eigen::Ref<Eigen::MatrixXd> MDiffXY) const;
1138
1143 bool CheckNLOVEV(const std::vector<double> &v) const;
1144
1148 virtual void Debugging(const std::vector<double> &input,
1149 std::vector<double> &output) const = 0;
1150
1157 const int &WhichMinimizer = Minimizer::WhichMinimizerDefault) const;
1158
1163 std::vector<std::vector<double>> SignSymmetries;
1168 void FindSignSymmetries();
1169
1173 void SetUseTreeLevel(bool val);
1174
1180 std::pair<std::vector<double>, std::vector<double>>
1181 initModel(std::string linestr);
1182
1191 std::vector<double> initModel(const std::vector<double> &par);
1192
1198 std::vector<double> resetScale(const double &newScale);
1199
1205 double CalculateRatioAlpha(const std::vector<double> &vev_symmetric,
1206 const std::vector<double> &vev_broken,
1207 const double &Temp) const;
1208
1213 Eigen::VectorXd NablaVCT(const std::vector<double> &v) const;
1218 Eigen::MatrixXd HessianCT(const std::vector<double> &v) const;
1219
1225 virtual std::vector<double> GetCTIdentities() const;
1226};
1227
1228} // namespace BSMPT
1229#endif /* CLASSPOTENTIALORIGIN_H_ */
The Class_Potential_Origin class Base class for all models. This class contains all numerical calcula...
Definition ClassPotentialOrigin.h:57
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:458
ModelID::ModelIDs get_Model() const
get_Model
Definition ClassPotentialOrigin.h:547
std::size_t nVEV
Definition ClassPotentialOrigin.h:142
const auto & Get_Curvature_Higgs_L2() const
Get_Curvature_Higgs_L2 allows const ref access to the L_{(S)}^{ij} Tensor.
Definition ClassPotentialOrigin.h:406
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:274
std::vector< std::vector< double > > HiggsRotationMatrix
Definition ClassPotentialOrigin.h:289
std::vector< std::vector< double > > SignSymmetries
Definition ClassPotentialOrigin.h:1163
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:476
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:212
double get_TripleHiggsCorrectionsCTPhysical(std::size_t i, std::size_t j, std::size_t k) const
get_TripleHiggsCorrectionsCTPhysical
Definition ClassPotentialOrigin.h:586
std::vector< double > get_parStored() const
get_parStored
Definition ClassPotentialOrigin.h:507
std::vector< double > vevTreeMin
Definition ClassPotentialOrigin.h:155
std::vector< double > get_parCTStored() const
get_parCTStored
Definition ClassPotentialOrigin.h:512
std::vector< std::vector< double > > DebyeGauge
DebyeGauge Stores the debye corrections to the mass matrix of the gauge bosons.
Definition ClassPotentialOrigin.h:384
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:295
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:853
std::vector< std::vector< std::vector< double > > > TripleHiggsCorrectionsCTPhysical
Definition ClassPotentialOrigin.h:178
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:220
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:329
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:491
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:323
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:78
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:311
std::vector< double > WeinbergForthDerivative() const
Definition ClassPotentialOrigin.cpp:1902
std::vector< std::complex< double > > QuarkMasses(const std::vector< double > &v) const
Definition ClassPotentialOrigin.cpp:3853
const ISMConstants SMConstants
SMConstants The SM constants used by the model.
Definition ClassPotentialOrigin.h:62
std::vector< std::vector< double > > Curvature_Higgs_CT_L2
Definition ClassPotentialOrigin.h:232
void CheckImplementation(const int &WhichMinimizer=Minimizer::WhichMinimizerDefault) const
Definition ClassPotentialOrigin.cpp:3544
std::vector< double > Curvature_Higgs_CT_L1
Definition ClassPotentialOrigin.h:228
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:151
virtual std::vector< double > calc_CT() const =0
std::vector< std::vector< std::vector< std::vector< double > > > > Curvature_Gauge_G2H2
Definition ClassPotentialOrigin.h:246
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:366
std::vector< std::vector< std::vector< double > > > TripleHiggsCorrectionsCWPhysical
Definition ClassPotentialOrigin.h:166
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:573
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:317
void CalculateDebye(bool forceCalculation=false)
Definition ClassPotentialOrigin.cpp:3223
std::vector< std::vector< std::vector< std::vector< double > > > > Curvature_Higgs_L4
Definition ClassPotentialOrigin.h:224
std::vector< std::vector< std::vector< double > > > Curvature_Higgs_CT_L3
Definition ClassPotentialOrigin.h:236
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:389
std::size_t nParCT
Definition ClassPotentialOrigin.h:82
void setUseIndexCol(std::string legend)
Definition ClassPotentialOrigin.cpp:3534
bool CalculatedTripleCopulings
CalculatedTripleCopulings Used to check if TripleHiggsCouplings has already been called.
Definition ClassPotentialOrigin.h:193
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:349
std::vector< std::vector< std::vector< std::complex< double > > > > Curvature_Lepton_F2H1
Definition ClassPotentialOrigin.h:260
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:3799
std::vector< std::vector< std::vector< double > > > TripleHiggsCorrectionsTreePhysical
Definition ClassPotentialOrigin.h:172
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:464
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:105
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:599
virtual void write() const =0
std::vector< std::vector< std::vector< std::complex< double > > > > Curvature_Quark_F2H1
Definition ClassPotentialOrigin.h:251
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:359
double scale
Definition ClassPotentialOrigin.h:73
std::vector< double > parStored
Definition ClassPotentialOrigin.h:87
bool CalcCouplingsdone
CalcCouplingsdone Used to check if CalculatePhysicalCouplings has already been called.
Definition ClassPotentialOrigin.h:188
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:130
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:3966
std::size_t get_NChargedHiggs() const
get_NChargedHiggs
Definition ClassPotentialOrigin.h:532
std::vector< std::vector< std::complex< double > > > Curvature_Lepton_F2
Definition ClassPotentialOrigin.h:264
ModelID::ModelIDs Model
Definition ClassPotentialOrigin.h:97
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:340
const std::vector< std::vector< double > > & get_DebyeHiggs() const
get_DebyeHiggs get the Debye corrections to the Higgs mass matrix
Definition ClassPotentialOrigin.h:553
std::size_t NNeutralHiggs
Definition ClassPotentialOrigin.h:101
void initVectors()
Definition ClassPotentialOrigin.cpp:3332
std::size_t get_NNeutralHiggs() const
get_NNeutralHiggs
Definition ClassPotentialOrigin.h:537
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:183
double FCW(double MassSquared) const
Definition ClassPotentialOrigin.cpp:80
int InputLineNumber
InputLineNumber Used for the error message in fbaseTri.
Definition ClassPotentialOrigin.h:198
Eigen::VectorXd NablaVCT(const std::vector< double > &v) const
NablaVCT.
Definition ClassPotentialOrigin.cpp:3986
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:300
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:284
std::vector< double > MassSquaredHiggs
MassSquaredHiggs Stores the masses of the Higgs Bosons calculated in CalculatePhysicalCouplings.
Definition ClassPotentialOrigin.h:269
double get_vevTreeMin(const std::size_t &k) const
get_vevTreeMin
Definition ClassPotentialOrigin.h:502
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:451
bool UseVCounterSimplified
UseVCounterSimplified Decides wether VCounterSimplified will be used or not. VCounterSimplified retur...
Definition ClassPotentialOrigin.h:864
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:434
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:444
bool UseTreeLevel
UseTreeLevel Enforces VEff to only use the tree-level potential.
Definition ClassPotentialOrigin.h:68
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:122
double get_scale() const
get_scale
Definition ClassPotentialOrigin.h:470
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:306
void SetUseTreeLevel(bool val)
Definition ClassPotentialOrigin.cpp:3746
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:379
std::vector< double > parCTStored
Definition ClassPotentialOrigin.h:92
std::size_t get_nPar() const
get_nPar
Definition ClassPotentialOrigin.h:481
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:344
const auto & Get_Curvature_Gauge_G2H2() const
Get_Curvature_Gauge_G2H2 allows const ref access to the G^{abij} Tensor.
Definition ClassPotentialOrigin.h:427
std::size_t NLepton
Definition ClassPotentialOrigin.h:137
virtual void set_gen(const std::vector< double > &par)=0
std::size_t NHiggs
Definition ClassPotentialOrigin.h:110
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:3927
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:241
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:354
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:335
std::vector< double > HiggsVev
Definition ClassPotentialOrigin.h:208
std::size_t get_nParCT() const
get_nParCT
Definition ClassPotentialOrigin.h:486
virtual std::vector< double > GetCTIdentities() const
GetCTIdentities.
Definition ClassPotentialOrigin.cpp:4032
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:279
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:3946
double fbaseTri(double MassSquaredA, double MassSquaredB, double MassSquaredC) const
Definition ClassPotentialOrigin.cpp:288
std::vector< double > get_vevTreeMin() const
get_vevTreeMin
Definition ClassPotentialOrigin.h:496
std::size_t get_NQuarks() const
get_NQuarks
Definition ClassPotentialOrigin.h:522
const auto & Get_Curvature_Higgs_L3() const
Get_Curvature_Higgs_L3 allows const ref access to the L_{(S)}^{ijk} Tensor.
Definition ClassPotentialOrigin.h:413
virtual bool CalculateDebyeGaugeSimplified()=0
std::size_t get_NHiggs() const
get_NHiggs
Definition ClassPotentialOrigin.h:527
std::vector< double > VEVSymmetric
Definition ClassPotentialOrigin.h:147
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:3752
Eigen::MatrixXd HessianCT(const std::vector< double > &v) const
HessianWeinberg.
Definition ClassPotentialOrigin.cpp:4010
std::vector< std::vector< double > > Curvature_Higgs_L2
Definition ClassPotentialOrigin.h:216
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:373
std::size_t NGauge
Definition ClassPotentialOrigin.h:116
const auto & Get_Curvature_Higgs_L4() const
Get_Curvature_Higgs_L4 allows const ref access to the L_{(S)}^{ijkl} Tensor.
Definition ClassPotentialOrigin.h:420
std::size_t get_NLepton() const
get_NLepton
Definition ClassPotentialOrigin.h:542
std::size_t get_NGauge() const
get_NGauge
Definition ClassPotentialOrigin.h:517
std::vector< std::vector< std::complex< double > > > Curvature_Quark_F2
Definition ClassPotentialOrigin.h:255
bool UseIndexCol
Definition ClassPotentialOrigin.h:203
std::vector< double > resetScale(const double &newScale)
resetScale changes the MSBar scale to newScale
Definition ClassPotentialOrigin.cpp:3787
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:3875
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:160
void FindSignSymmetries()
FindSignSymmetries checks for all possible sign changes in the VEV components and checks for all poss...
Definition ClassPotentialOrigin.cpp:3710
void set_InputLineNumber(int InputLineNumber_in)
set_InputLineNumber
Definition ClassPotentialOrigin.h:562
ModelIDs
The ModelIDs enum containing all IDs for identifying the Models.
Definition IncludeAllModels.h:32
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:32
const double C_CWcbHiggs
C_CWcbHiggs constant used in the CW potential for Higgs bosons.
Definition ClassPotentialOrigin.h:49
const double C_CWcbGB
C_CWcbGB constant used in the CW potential for gauge bosons.
Definition ClassPotentialOrigin.h:45
const double C_threshold
C_threshold threshold to check if a mass is numerically zero.
Definition ClassPotentialOrigin.h:37
const bool C_UseParwani
C_UseParwani Use the Parwani Method instead of Arnold-Espinosa.
Definition ClassPotentialOrigin.h:27
The ISMConstants struct containing all necessary SM constants.
Definition SMparam.h:24
Definition transition_tracer.h:143