BSMPT 3.0.7
BSMPT - Beyond the Standard Model Phase Transitions : A C++ package for the computation of the EWPT in BSM models
|
The Class_Potential_Origin class Base class for all models. This class contains all numerical calculations on the tensors and the inherited classes only have to set them. More...
#include <BSMPT/models/ClassPotentialOrigin.h>
Public Member Functions | |
Class_Potential_Origin (const ISMConstants &smConstants) | |
const auto & | Get_Curvature_Higgs_L2 () const |
Get_Curvature_Higgs_L2 allows const ref access to the L_{(S)}^{ij} Tensor. | |
const auto & | Get_Curvature_Higgs_L3 () const |
Get_Curvature_Higgs_L3 allows const ref access to the L_{(S)}^{ijk} Tensor. | |
const auto & | Get_Curvature_Higgs_L4 () const |
Get_Curvature_Higgs_L4 allows const ref access to the L_{(S)}^{ijkl} Tensor. | |
const auto & | Get_Curvature_Gauge_G2H2 () const |
Get_Curvature_Gauge_G2H2 allows const ref access to the G^{abij} Tensor. | |
const auto & | Get_Curvature_Lepton_F2H1 () const |
Get_Curvature_Lepton_F2H1 allows const ref access to the Y^{IJk} Tensor for Leptons. | |
const auto & | Get_Curvature_Lepton_F2 () const |
Get_Curvature_Lepton_F2 allows const ref access to the Y^{IJ} Tensor for Leptons. | |
const auto & | Get_Curvature_Quark_F2H1 () const |
Get_Curvature_Quark_F2H1 allows const ref access to the Y^{IJk} Tensor for Quarks. | |
const auto & | Get_Curvature_Quark_F2 () const |
Get_Curvature_Quark_F2 allows const ref access to the Y^{IJ} Tensor for Quarks. | |
const std::vector< std::size_t > & | Get_VevOrder () const |
Get_VevOrder allows const ref access to the VevOrder vector. | |
double | get_scale () const |
get_scale | |
void | set_scale (double scale_new) |
set_scale sets the MSBar renormalisation scale to scale_new | |
std::size_t | get_nPar () const |
get_nPar | |
std::size_t | get_nParCT () const |
get_nParCT | |
std::size_t | get_nVEV () const |
get_nVEV | |
std::vector< double > | get_vevTreeMin () const |
get_vevTreeMin | |
double | get_vevTreeMin (const std::size_t &k) const |
get_vevTreeMin | |
std::vector< double > | get_parStored () const |
get_parStored | |
std::vector< double > | get_parCTStored () const |
get_parCTStored | |
std::size_t | get_NGauge () const |
get_NGauge | |
std::size_t | get_NQuarks () const |
get_NQuarks | |
std::size_t | get_NHiggs () const |
get_NHiggs | |
std::size_t | get_NChargedHiggs () const |
get_NChargedHiggs | |
std::size_t | get_NNeutralHiggs () const |
get_NNeutralHiggs | |
std::size_t | get_NLepton () const |
get_NLepton | |
ModelID::ModelIDs | get_Model () const |
get_Model | |
const std::vector< std::vector< double > > & | get_DebyeHiggs () const |
get_DebyeHiggs get the Debye corrections to the Higgs mass matrix | |
void | set_InputLineNumber (int InputLineNumber_in) |
set_InputLineNumber | |
double | get_TripleHiggsCorrectionsTreePhysical (std::size_t i, std::size_t j, std::size_t k) const |
get_TripleHiggsCorrectionsTreePhysical | |
double | get_TripleHiggsCorrectionsCTPhysical (std::size_t i, std::size_t j, std::size_t k) const |
get_TripleHiggsCorrectionsCTPhysical | |
double | get_TripleHiggsCorrectionsCWPhysical (std::size_t i, std::size_t j, std::size_t k) const |
get_TripleHiggsCorrectionsCWPhysical | |
void | setUseIndexCol (std::string legend) |
bool | getUseIndexCol () |
void | sym2Dim (std::vector< std::vector< double > > &Tensor2Dim, std::size_t Nk1, std::size_t Nk2) |
sym2Dim Symmetrize scalar 2-dim tensor | |
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 | |
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 | |
void | initVectors () |
virtual double | VEff (const std::vector< double > &v, double Temp=0, int diff=0, int Order=1) const |
double | VTree (const std::vector< double > &v, int diff=0, bool ForceExplicitCalculation=false) const |
double | CounterTerm (const std::vector< double > &v, int diff=0, bool ForceExplicitCalculation=false) const |
double | V1Loop (const std::vector< double > &v, double Temp, int diff) const |
virtual double | EWSBVEV (const std::vector< double > &v) const |
virtual void | SetEWVEVZero (std::vector< double > &sol) const |
SetEWVEVZero Set all VEV directions in sol-vector to zero that contibute to EW VEV. | |
virtual void | ReadAndSet (const std::string &linestr, std::vector< double > &par)=0 |
virtual std::vector< std::string > | addLegendCT () const =0 |
virtual std::vector< std::string > | addLegendTemp () const =0 |
virtual std::vector< std::string > | addLegendTripleCouplings () const =0 |
virtual std::vector< std::string > | addLegendVEV () const =0 |
virtual void | set_gen (const std::vector< double > &par)=0 |
virtual void | set_CT_Pot_Par (const std::vector< double > &par)=0 |
virtual void | write () const =0 |
void | set_All (const std::vector< double > &par, const std::vector< double > &parCT) |
virtual void | SetCurvatureArrays ()=0 |
void | CalculatePhysicalCouplings () |
std::vector< double > | WeinbergFirstDerivative () const |
std::vector< double > | WeinbergSecondDerivative () const |
Eigen::MatrixXd | WeinbergSecondDerivativeAsMatrixXd () const |
std::vector< double > | WeinbergThirdDerivative () const |
std::vector< double > | WeinbergForthDerivative () const |
void | CalculateDebye (bool forceCalculation=false) |
void | CalculateDebyeGauge () |
void | Prepare_Triple () |
virtual bool | CalculateDebyeSimplified ()=0 |
virtual bool | CalculateDebyeGaugeSimplified ()=0 |
virtual double | VTreeSimplified (const std::vector< double > &v) const =0 |
virtual double | VCounterSimplified (const std::vector< double > &v) const =0 |
std::vector< double > | HiggsMassesSquared (const std::vector< double > &v, const double &Temp=0, const int &diff=0) const |
Eigen::MatrixXd | HiggsMassMatrix (const std::vector< double > &v, double Temp=0, int diff=0) const |
HiggsMassMatrix calculates the Higgs mass matrix. | |
std::vector< double > | GaugeMassesSquared (const std::vector< double > &v, const double &Temp=0, const int &diff=0) const |
std::vector< double > | QuarkMassesSquared (const std::vector< double > &v, const int &diff=0) const |
std::vector< double > | LeptonMassesSquared (const std::vector< double > &v, const int &diff=0) const |
std::vector< std::complex< double > > | QuarkMasses (const std::vector< double > &v) const |
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 $. | |
std::vector< std::complex< double > > | LeptonMasses (const std::vector< double > &v) const |
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_k $. | |
virtual void | TripleHiggsCouplings ()=0 |
double | fbase (double MassSquaredA, double MassSquaredB) const |
double | fbaseTri (double MassSquaredA, double MassSquaredB, double MassSquaredC) const |
double | fbaseFour (double MassSquaredA, double MassSquaredB, double MassSquaredC, double MassSquaredD) const |
virtual std::vector< double > | calc_CT () const =0 |
double | CWTerm (double MassSquared, double cb, int diff=0) const |
double | FCW (double MassSquared) const |
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.o.f. std::size_to account. | |
double | boson_legacy (double MassSquared, double Temp, double cb, int diff=0) const |
double | fermion (double MassSquared, double Temp, int diff=0) const |
Calculation of the fermionic std::size_tegral + Coleman-Weinberg potential without taking d.o.f. std::size_to account. | |
double | fermion_legacy (double Mass, double Temp, int diff=0) const |
double | Vl (double MassSquared, double Temp, int n, int diff=0) const |
double | Vsf (double MassSquared, double Temp, int n, int diff=0) const |
double | Vsb (double MassSquaredIn, double Temp, int n, int diff=0) const |
void | resetbools () |
std::vector< double > | MinimizeOrderVEV (const std::vector< double > &VEVminimizer) const |
std::vector< double > | FirstDerivativeOfEigenvalues (const Eigen::Ref< Eigen::MatrixXcd > M, const Eigen::Ref< Eigen::MatrixXcd > MDiff) const |
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 |
bool | CheckNLOVEV (const std::vector< double > &v) const |
virtual void | Debugging (const std::vector< double > &input, std::vector< double > &output) const =0 |
void | CheckImplementation (const int &WhichMinimizer=Minimizer::WhichMinimizerDefault) const |
void | FindSignSymmetries () |
FindSignSymmetries checks for all possible sign changes in the VEV components and checks for all possible Z2 symmetries. | |
void | SetUseTreeLevel (bool val) |
std::pair< std::vector< double >, std::vector< double > > | initModel (std::string linestr) |
std::vector< double > | initModel (const std::vector< double > &par) |
std::vector< double > | resetScale (const double &newScale) |
resetScale changes the MSBar scale to newScale | |
double | CalculateRatioAlpha (const std::vector< double > &vev_symmetric, const std::vector< double > &vev_broken, const double &Temp) const |
Eigen::VectorXd | NablaVCT (const std::vector< double > &v) const |
NablaVCT. | |
Eigen::MatrixXd | HessianCT (const std::vector< double > &v) const |
HessianWeinberg. | |
virtual std::vector< double > | GetCTIdentities () const |
GetCTIdentities. | |
Public Attributes | |
const ISMConstants | SMConstants |
SMConstants The SM constants used by the model. | |
bool | UseVTreeSimplified = false |
UseVTreeSimplified Decides wether VTreeSimplified will be used or not. VTreeSimplified returns 0 if UseVTreeSimplified is false Set in constructor of the implemented models. | |
bool | UseVCounterSimplified = false |
UseVCounterSimplified Decides wether VCounterSimplified will be used or not. VCounterSimplified returns 0 if UseVCounterSimplified is false Set in constructor of the implemented models. | |
std::vector< std::vector< double > > | SignSymmetries |
Protected Attributes | |
bool | UseTreeLevel = false |
UseTreeLevel Enforces VEff to only use the tree-level potential. | |
double | scale |
std::size_t | nPar = 0 |
std::size_t | nParCT = 0 |
std::vector< double > | parStored |
std::vector< double > | parCTStored |
ModelID::ModelIDs | Model = ModelID::ModelIDs::NotSet |
std::size_t | NNeutralHiggs = 0 |
std::size_t | NChargedHiggs = 0 |
std::size_t | NHiggs = NNeutralHiggs + NChargedHiggs |
std::size_t | NGauge = 4 |
std::size_t | NQuarks = 12 |
std::size_t | NColour = 3 |
NColour Number of colours of the quarks. Do not change this in the current version as we only investigate extended Higgs sectors. If you want to extend the other sectors as well the Debye corrections have to be calculated by hand. | |
std::size_t | NLepton = 9 |
std::size_t | nVEV = 0 |
std::vector< double > | VEVSymmetric |
std::vector< double > | vevTree |
std::vector< double > | vevTreeMin |
std::vector< double > | TripleHiggsCorrectionsCW |
std::vector< std::vector< std::vector< double > > > | TripleHiggsCorrectionsCWPhysical |
std::vector< std::vector< std::vector< double > > > | TripleHiggsCorrectionsTreePhysical |
std::vector< std::vector< std::vector< double > > > | TripleHiggsCorrectionsCTPhysical |
bool | SetCurvatureDone = false |
SetCurvatureDone Used to check if the tensors are set. | |
bool | CalcCouplingsdone = false |
CalcCouplingsdone Used to check if CalculatePhysicalCouplings has already been called. | |
bool | CalculatedTripleCopulings = false |
CalculatedTripleCopulings Used to check if TripleHiggsCouplings has already been called. | |
int | InputLineNumber = -1 |
InputLineNumber Used for the error message in fbaseTri. | |
bool | UseIndexCol = false |
std::vector< double > | HiggsVev |
std::vector< double > | Curvature_Higgs_L1 |
std::vector< std::vector< double > > | Curvature_Higgs_L2 |
std::vector< std::vector< std::vector< double > > > | Curvature_Higgs_L3 |
std::vector< std::vector< std::vector< std::vector< double > > > > | Curvature_Higgs_L4 |
std::vector< double > | Curvature_Higgs_CT_L1 |
std::vector< std::vector< double > > | Curvature_Higgs_CT_L2 |
std::vector< std::vector< std::vector< double > > > | Curvature_Higgs_CT_L3 |
std::vector< std::vector< std::vector< std::vector< double > > > > | Curvature_Higgs_CT_L4 |
std::vector< std::vector< std::vector< std::vector< double > > > > | Curvature_Gauge_G2H2 |
std::vector< std::vector< std::vector< std::complex< double > > > > | Curvature_Quark_F2H1 |
std::vector< std::vector< std::complex< double > > > | Curvature_Quark_F2 |
std::vector< std::vector< std::vector< std::complex< double > > > > | Curvature_Lepton_F2H1 |
std::vector< std::vector< std::complex< double > > > | Curvature_Lepton_F2 |
std::vector< double > | MassSquaredHiggs |
MassSquaredHiggs Stores the masses of the Higgs Bosons calculated in CalculatePhysicalCouplings. | |
std::vector< double > | MassSquaredGauge |
MassSquaredGauge Stores the masses of the gauge Bosons calculated in CalculatePhysicalCouplings. | |
std::vector< double > | MassSquaredQuark |
MassSquaredQuark Stores the masses of the quarks calculated in CalculatePhysicalCouplings. | |
std::vector< double > | MassSquaredLepton |
MassSquaredLepton Stores the masses of the leptons calculated in CalculatePhysicalCouplings. | |
std::vector< std::vector< double > > | HiggsRotationMatrix |
std::vector< std::vector< std::vector< std::vector< double > > > > | Couplings_Higgs_Quartic |
Couplings_Higgs_Quartic Stores the quartic Higgs couplings in the mass base. | |
std::vector< std::vector< std::vector< double > > > | Couplings_Higgs_Triple |
Couplings_Higgs_Triple Stores the triple Higgs couplings in the mass base. | |
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. | |
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. | |
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 base. | |
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. | |
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 base. | |
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 base. | |
std::vector< std::vector< std::vector< double > > > | LambdaGauge_3 |
LambdaGauge_3 Stores the Lambda_{(G)}^{abi} tensor. | |
std::vector< std::vector< std::vector< double > > > | LambdaHiggs_3 |
LambdaHiggs_3 Stores the Lambda_{(S)}^{ijk} tensor. | |
std::vector< std::vector< std::vector< double > > > | LambdaHiggs_3_CT |
LambdaHiggs_3 Stores the Lambda_{(S)}^{ijk} tensor for the counterterm parameters. | |
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_{(F)}^{IJ} w.r.t. the Higgs field k. | |
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_{(F)}^{IJ} w.r.t. the Higgs field k. | |
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_{(F)}^{IJ} w.r.t. the Higgs fields k and m. | |
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 Lambda_{(F)}^{IJ} w.r.t. the Higgs fields k and m. | |
std::vector< std::vector< double > > | DebyeHiggs |
DebyeHiggs Stores the debye corrections to the mass matrix of the Higgs bosons. | |
std::vector< std::vector< double > > | DebyeGauge |
DebyeGauge Stores the debye corrections to the mass matrix of the gauge bosons. | |
std::vector< std::size_t > | VevOrder |
VevOrder Stores the matching order used in MinimizeOrderVEV, set in the constructor of the model. | |
The Class_Potential_Origin class Base class for all models. This class contains all numerical calculations on the tensors and the inherited classes only have to set them.
|
pure virtual |
Adds the name of the counterterms to the legend of the output file. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
pure virtual |
Adds the name of the VEVs, together with the critical VEV, the critical temperature and the ratio, to the legend of the output file. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
pure virtual |
Adds the name of the triple Higgs couplings at T=0 to the legend of the output file. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
pure virtual |
Adds the name of the VEVs to the legend of the output file. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
double BSMPT::Class_Potential_Origin::boson | ( | double | MassSquared, |
double | Temp, | ||
double | cb, | ||
int | diff = 0 |
||
) | const |
Calculation of the bosonic std::size_tegral + Coleman-Weinberg potential without taking d.o.f. std::size_to account.
MassSquared | = m^2 of the particle |
Temp | = Temperature at which the Debye masses and std::size_tegrals should be calculated |
cb | = Parameter of the renormalisation-Scheme in the Coleman-Weinberg potential |
diff | 0 returns the value of the std::size_tegral and diff >0 the derivative w.r.t. m^2 and diff = -1 w.r.t. Temp |
double BSMPT::Class_Potential_Origin::boson_legacy | ( | double | MassSquared, |
double | Temp, | ||
double | cb, | ||
int | diff = 0 |
||
) | const |
Deprecated version of boson() as present in the v1.X release. Still here for legacy reasons
|
pure virtual |
Calculates the counterterm parameters. Here you need to work out the scheme and implement the formulas. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
void BSMPT::Class_Potential_Origin::CalculateDebye | ( | bool | forceCalculation = false | ) |
Calculates the Debye corrections to the Higgs mass matrix. If you can provide CalculateDebyeSimplified() with the Matrix as this will reduce the runtime.
forceCalculation | Forces the caclulation, even if the model implements a simplified version. The simplified calculation will be skipped. |
void BSMPT::Class_Potential_Origin::CalculateDebyeGauge | ( | ) |
Calculates the Debye corrections to the gauge sector. By using CalculateDebyeGaugeSimplified() the runtime can be reduced.
|
pure virtual |
You can give the explicit Debye corrections to the gauge boson mass matrix with this function if you know it. Otherwise it is calculated via CalculateDebyeGauge(). If you know the corrections use this and let the function return true, this will save you a lot of computing time.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
pure virtual |
You can give the explicit Debye corrections to the Higgs mass matrix with this function if you know it. Otherwise it is calculated via CalculateDebye(). If you know the corrections use this and let the function return true, this will save you a lot of computing time.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
void BSMPT::Class_Potential_Origin::CalculatePhysicalCouplings | ( | ) |
Calculates all triple and quartic couplings in the physical basis
double BSMPT::Class_Potential_Origin::CalculateRatioAlpha | ( | const std::vector< double > & | vev_symmetric, |
const std::vector< double > & | vev_broken, | ||
const double & | Temp | ||
) | const |
Calculate the ratio of the latent heat w.r.t. the energy density of the plasma background Formula taken from
void BSMPT::Class_Potential_Origin::CheckImplementation | ( | const int & | WhichMinimizer = Minimizer::WhichMinimizerDefault | ) | const |
Checks if the tensors are correctly implemented. For this the fermion, quark and gauge boson masses are calculated and printed next to the values defined in SMparah.h
bool BSMPT::Class_Potential_Origin::CheckNLOVEV | ( | const std::vector< double > & | v | ) | const |
This function will check if the VEV at NLO is still close enough to the LO VEV. v has to be given in the configuration of the minimizer.
double BSMPT::Class_Potential_Origin::CounterTerm | ( | const std::vector< double > & | v, |
int | diff = 0 , |
||
bool | ForceExplicitCalculation = false |
||
) | const |
Calculates the counterterm potential and its derivatives
v | the configuration of all VEVs at which the potential should be calculated |
diff | 0 returns the potential and i!= 0 returns the derivative of the potential w.r.t v_i |
ForceExplicitCalculation | Calculate the tensors directly from the tensors even if VCounterSimplified() is given |
double BSMPT::Class_Potential_Origin::CWTerm | ( | double | MassSquared, |
double | cb, | ||
int | diff = 0 |
||
) | const |
Calculates the Coleman-Weinberg contribution of a particle with m^2 = MassSquared and the constant scheme-dependent parameter cb as well as the first derivative with respect to m^2.
|
pure virtual |
This is a possible debugging function.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
virtual |
This function calculates the EW breaking VEV from all contributing field configurations.
double BSMPT::Class_Potential_Origin::fbase | ( | double | MassSquaredA, |
double | MassSquaredB | ||
) | const |
Calculates the function f_{ab}^{(1)} (see https://arxiv.org/abs/1606.07069 eq 3.9) needed for the derivatives of the Coleman Weinberg potential.
double BSMPT::Class_Potential_Origin::fbaseFour | ( | double | MassSquaredA, |
double | MassSquaredB, | ||
double | MassSquaredC, | ||
double | MassSquaredD | ||
) | const |
Calculates the function f_{abcd}^{(1)} (see https://arxiv.org/abs/1606.07069 eq 3.9) needed for the 4th derivatives of the Coleman Weinberg potential.
double BSMPT::Class_Potential_Origin::fbaseTri | ( | double | MassSquaredA, |
double | MassSquaredB, | ||
double | MassSquaredC | ||
) | const |
Calculates the function f_{abc}^{(1)} (see https://arxiv.org/abs/1606.07069 eq 3.9) needed for the 3rd derivatives of the Coleman Weinberg potential.
double BSMPT::Class_Potential_Origin::FCW | ( | double | MassSquared | ) | const |
Calculates Re(log(MassSquared)) and returns 0 if the argument is too small as this function is only called with an (m^2)^n in front of it.
double BSMPT::Class_Potential_Origin::fermion | ( | double | MassSquared, |
double | Temp, | ||
int | diff = 0 |
||
) | const |
Calculation of the fermionic std::size_tegral + Coleman-Weinberg potential without taking d.o.f. std::size_to account.
MassSquared | = m^2 of the particle |
Temp | = temperature at which the std::size_tegrals should be calculated |
diff | : 0 = Value of the potential, diff > 0 returns the derivative w.r.t. m^2 and diff=-1 w.r.t Temp |
double BSMPT::Class_Potential_Origin::fermion_legacy | ( | double | Mass, |
double | Temp, | ||
int | diff = 0 |
||
) | const |
Deprecated version of fermion() as present in the v1.X release. Still here for legacy reasons
std::vector< double > BSMPT::Class_Potential_Origin::FirstDerivativeOfEigenvalues | ( | const Eigen::Ref< Eigen::MatrixXcd > | M, |
const Eigen::Ref< Eigen::MatrixXcd > | MDiff | ||
) | const |
Calculates the first derivatives of the eigenvalues of a given matrix
M | : the original matrix |
MDiff | : the element-wise first derivative of the matrix M with respect to the parameter you want to consider |
std::vector< double > BSMPT::Class_Potential_Origin::GaugeMassesSquared | ( | const std::vector< double > & | v, |
const double & | Temp = 0 , |
||
const int & | diff = 0 |
||
) | const |
Calculates the gauge mass matrix and saves all eigenvalues
v | the configuration of all VEVs at which the eigenvalues should be evaluated |
Temp | The temperature at which the Debye corrected masses should be calculated |
diff | 0 returns the masses and i!=0 returns the derivative of m^2 w.r.t v_i, -1 returns the derivative w.r.t. the temperature |
|
inline |
Get_Curvature_Gauge_G2H2 allows const ref access to the G^{abij} Tensor.
|
inline |
Get_Curvature_Higgs_L2 allows const ref access to the L_{(S)}^{ij} Tensor.
|
inline |
Get_Curvature_Higgs_L3 allows const ref access to the L_{(S)}^{ijk} Tensor.
|
inline |
Get_Curvature_Higgs_L4 allows const ref access to the L_{(S)}^{ijkl} Tensor.
|
inline |
Get_Curvature_Lepton_F2 allows const ref access to the Y^{IJ} Tensor for Leptons.
|
inline |
Get_Curvature_Lepton_F2H1 allows const ref access to the Y^{IJk} Tensor for Leptons.
|
inline |
Get_Curvature_Quark_F2 allows const ref access to the Y^{IJ} Tensor for Quarks.
|
inline |
Get_Curvature_Quark_F2H1 allows const ref access to the Y^{IJk} Tensor for Quarks.
|
inline |
get_DebyeHiggs get the Debye corrections to the Higgs mass matrix
|
inline |
get_Model
|
inline |
get_NChargedHiggs
|
inline |
get_NGauge
|
inline |
get_NHiggs
|
inline |
get_NLepton
|
inline |
get_NNeutralHiggs
|
inline |
get_nPar
|
inline |
get_nParCT
|
inline |
get_NQuarks
|
inline |
get_nVEV
|
inline |
get_parCTStored
|
inline |
get_parStored
|
inline |
get_scale
|
inline |
get_TripleHiggsCorrectionsCTPhysical
i | |
j | |
k |
|
inline |
get_TripleHiggsCorrectionsCWPhysical
i | |
j | |
k |
|
inline |
get_TripleHiggsCorrectionsTreePhysical
i | |
j | |
k |
|
inline |
Get_VevOrder allows const ref access to the VevOrder vector.
|
inline |
get_vevTreeMin
|
inline |
get_vevTreeMin
k |
|
virtual |
GetCTIdentities.
bool BSMPT::Class_Potential_Origin::getUseIndexCol | ( | ) |
get UseIndexCol var
Eigen::MatrixXd BSMPT::Class_Potential_Origin::HessianCT | ( | const std::vector< double > & | v | ) | const |
HessianWeinberg.
std::vector< double > BSMPT::Class_Potential_Origin::HiggsMassesSquared | ( | const std::vector< double > & | v, |
const double & | Temp = 0 , |
||
const int & | diff = 0 |
||
) | const |
Calculates the Higgs mass matrix and saves all eigenvalues
v | the configuration of all VEVs at which the eigenvalues should be evaluated |
Temp | The temperature at which the Debye corrected masses should be calculated |
diff | 0 returns the masses and i!=0 returns the derivative of m^2 w.r.t v_i, i = -1 returns the derivative w.r.t. to the temperature |
MatrixXd BSMPT::Class_Potential_Origin::HiggsMassMatrix | ( | const std::vector< double > & | v, |
double | Temp = 0 , |
||
int | diff = 0 |
||
) | const |
HiggsMassMatrix calculates the Higgs mass matrix.
v | the configuration of all VEVs at which the Mass Matrix should be evaluated |
Temp | The temperature at which the Debye corrected masses should be calculated |
diff | 0 returns the masses and i!=0 returns the derivative the Mass Matrix w.r.t v_i, i = -1 returns the derivative w.r.t. to the temperature |
std::vector< double > BSMPT::Class_Potential_Origin::initModel | ( | const std::vector< double > & | par | ) |
Gets the parameter in a vector in the same way as in set_gen and calls resetbools, set_gen, calc_CT, set_CT_Pot_Par,CalculateDebye and CalculateDebyeGauge
par | Parameters to define the parameter point. Same input as in set_gen() |
std::pair< std::vector< double >, std::vector< double > > BSMPT::Class_Potential_Origin::initModel | ( | std::string | linestr | ) |
Gets the parameter line as an Input and calls resetbools, ReadAndSet, calc_CT, set_CT_Pot_Par,CalculateDebye and CalculateDebyeGauge
void BSMPT::Class_Potential_Origin::initVectors | ( | ) |
Initializes all vectors needed for the calculations.
std::vector< std::complex< double > > BSMPT::Class_Potential_Origin::LeptonMasses | ( | const std::vector< double > & | v | ) | const |
Calculates the quark mass matrix and saves all eigenvalues, this assumes the same masses for different colours.
v | the configuration of all VEVs at which the eigenvalues should be evaluated |
std::vector< double > BSMPT::Class_Potential_Origin::LeptonMassesSquared | ( | const std::vector< double > & | v, |
const int & | diff = 0 |
||
) | const |
Calculates the lepton mass matrix and saves all eigenvalues
v | the configuration of all VEVs at which the eigenvalues should be evaluated |
diff | 0 returns the masses and i!=0 returns the derivative of m^2 w.r.t v_i |
MatrixXcd BSMPT::Class_Potential_Origin::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_k $.
v | the configuration of all VEVs at which the matrix should be calculated |
std::vector< double > BSMPT::Class_Potential_Origin::MinimizeOrderVEV | ( | const std::vector< double > & | VEVminimizer | ) | const |
This will convert the vector VEVminimizer with the VEVs used by the minimizer to a vector with the complete VEV configuration VEVFunction used by the Class functions. The mapping is done with the VevOrder vector which has to be set in the constructor of your model
VEVminimizer | vector with the VEVs, used as inputs for the minimizers |
Eigen::VectorXd BSMPT::Class_Potential_Origin::NablaVCT | ( | const std::vector< double > & | v | ) | const |
NablaVCT.
void BSMPT::Class_Potential_Origin::Prepare_Triple | ( | ) |
Sets a tensor needed to calculate the contribution of the counterterm potential to the triple Higgs couplings.
std::vector< std::complex< double > > BSMPT::Class_Potential_Origin::QuarkMasses | ( | const std::vector< double > & | v | ) | const |
Calculates the quark mass matrix and saves all eigenvalues, this assumes the same masses for different colours.
v | the configuration of all VEVs at which the eigenvalues should be evaluated |
std::vector< double > BSMPT::Class_Potential_Origin::QuarkMassesSquared | ( | const std::vector< double > & | v, |
const int & | diff = 0 |
||
) | const |
Calculates the quark mass matrix and saves all eigenvalues, this assumes the same masses for different colours.
v | the configuration of all VEVs at which the eigenvalues should be evaluated |
diff | 0 returns the masses and i!=0 returns the derivative of m^2 w.r.t v_i |
Eigen::MatrixXcd BSMPT::Class_Potential_Origin::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 $.
v | the configuration of all VEVs at which the matrix should be calculated |
|
pure virtual |
Reads the string linestr and sets the parameter point
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
void BSMPT::Class_Potential_Origin::resetbools | ( | ) |
Resets all bools. Needed if you want to deal with multiple posize_ts one after another with the same pointer.
std::vector< double > BSMPT::Class_Potential_Origin::resetScale | ( | const double & | newScale | ) |
resetScale changes the MSBar scale to newScale
newScale | the new Scale in GeV |
std::vector< double > BSMPT::Class_Potential_Origin::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 |
This function calculates the second derivatives of all eigenvalues. The matrix must not have a repeated eigenvalue for this!
M | : the original matrix |
MDiffX | : the element-wise first derivative of the matrix M with respect to the first parameter you want to consider |
MDiffY | : the element-wise first derivative of the matrix M with respect to the second parameter you want to consider |
MDiffXY | : the element-wise second derivative of the matrix M with respect to both parameters you want to consider |
void BSMPT::Class_Potential_Origin::set_All | ( | const std::vector< double > & | par, |
const std::vector< double > & | parCT | ||
) |
This will call set_gen(par), SetCurvatureArrays, set_CT_Pot_Par(parCT), CalculateDebye() as well as CalculateDebyeGauge()
|
pure virtual |
Reads the counterterm parameters from the vector 'par' and sets them to the model parameters as well as the Tensors Curvature_Higgs_CT_L1
, Curvature_Higgs_CT_L2
, Curvature_Higgs_CT_L3
, Curvature_Higgs_CT_L4
. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
pure virtual |
Reads the Lagrangian parameters from the vector 'par' and sets them to the model parameters. This also sets the VEV configuration as well as the scale. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
inline |
set_InputLineNumber
InputLineNumber_in | value to set InputLineNumber |
|
inline |
set_scale sets the MSBar renormalisation scale to scale_new
scale_new |
|
pure virtual |
This will set all the tensors needed to describe the tree-level Lagrangian except for the counterterms in the potential. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
virtual |
SetEWVEVZero Set all VEV directions in sol-vector to zero that contibute to EW VEV.
sol | solution in nVEV-space |
void BSMPT::Class_Potential_Origin::setUseIndexCol | ( | std::string | legend | ) |
Sets the UseIndexCol var
void BSMPT::Class_Potential_Origin::SetUseTreeLevel | ( | bool | val | ) |
Set the parameter UseTreeLevel to the input
void BSMPT::Class_Potential_Origin::sym2Dim | ( | std::vector< std::vector< double > > & | Tensor2Dim, |
std::size_t | Nk1, | ||
std::size_t | Nk2 | ||
) |
sym2Dim Symmetrize scalar 2-dim tensor
Tensor2Dim | 2-dim scalar tensor |
Nk1 | number of first indices |
Nk2 | number of second indices |
void BSMPT::Class_Potential_Origin::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
Tensor3Dim | 4-dim scalar tensor |
Nk1 | number of first indices |
Nk2 | number of second indices |
Nk3 | number of third indices |
void BSMPT::Class_Potential_Origin::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
Tensor4Dim | 4-dim scalar tensor |
Nk1 | number of first indices |
Nk2 | number of second indices |
Nk3 | number of third indices |
Nk4 | number of forth indices |
|
pure virtual |
Calculates the triple Higgs couplings at NLO in the mass basis.
Use the vector TripleHiggsCorrectionsCWPhysical to save your couplings and set the nTripleCouplings to the number of couplings you want as output.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
double BSMPT::Class_Potential_Origin::V1Loop | ( | const std::vector< double > & | v, |
double | Temp, | ||
int | diff | ||
) | const |
Calculates the Coleman-Weinberg and temperature-dependent 1-loop part of the effective potential and its derivatives.
diff | 0: No derivative, > 0 derivative w.r.t to the Higgs field; -1 : derivative w.r.t to the temperature |
v | the configuration of all VEVs at which the potential should be calculated |
Temp | the temperature at which the potential should be evaluated |
In case of diff != 0 the mass vectors will directly have the derivatives of the masses
|
pure virtual |
You can give the explicit form of your counterterm potential here. This speeds up the computation time.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
virtual |
Calculates the effective potential and its derivatives.
v | vev configuration at which the potential should be evaluated |
Temp | temperature at which the potential should be evaluated |
diff | Switch for the derivative of the potential. Default is 0 for the value of the potential |
Order | 0 returns the tree level potential and 1 the NLO potential. Default value is the NLO potential |
double BSMPT::Class_Potential_Origin::Vl | ( | double | MassSquared, |
double | Temp, | ||
int | n, | ||
int | diff = 0 |
||
) | const |
Calculates the large m^2/T^2 approximation to order n for the temperature-dependent std::size_tegrals.
double BSMPT::Class_Potential_Origin::Vsb | ( | double | MassSquaredIn, |
double | Temp, | ||
int | n, | ||
int | diff = 0 |
||
) | const |
Calculates the small m^2/T^2 approximation of the bosonic temperature dependent std::size_tegral to the n-th order.
double BSMPT::Class_Potential_Origin::Vsf | ( | double | MassSquared, |
double | Temp, | ||
int | n, | ||
int | diff = 0 |
||
) | const |
Calculates the small m^2/T^2 approximation of the fermionic temperature dependent std::size_tegral to the n-th order.
MassSquared | m^2 |
Temp | temperature |
n | the order to which expand the expansion |
diff | 0 = value; >0 derivative w.r.t. |
double BSMPT::Class_Potential_Origin::VTree | ( | const std::vector< double > & | v, |
int | diff = 0 , |
||
bool | ForceExplicitCalculation = false |
||
) | const |
Calculates the tree-level potential and its derivatives.
v | the configuration of all VEVs at which the potential should be calculated |
diff | 0 returns the potential and i!= 0 returns the derivative of the potential w.r.t v_i |
ForceExplicitCalculation | Calculate the tensors directly from the tensors even if VTreeSimplified() is given |
|
pure virtual |
You can give the explicit form of your tree-level potential here. This speeds up the computation time.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
std::vector< double > BSMPT::Class_Potential_Origin::WeinbergFirstDerivative | ( | ) | const |
Calculates the first derivative of the Coleman-Weinberg potential evaluated at the tree-level minimum.
std::vector< double > BSMPT::Class_Potential_Origin::WeinbergForthDerivative | ( | ) | const |
Calculates the forth derivative of the Coleman-Weinberg potential at the tree-level minimum.
std::vector< double > BSMPT::Class_Potential_Origin::WeinbergSecondDerivative | ( | ) | const |
Calculates the second derivative of the Coleman-Weinberg potential at the tree-level minimum.
Eigen::MatrixXd BSMPT::Class_Potential_Origin::WeinbergSecondDerivativeAsMatrixXd | ( | ) | const |
Calculates the second derivative of the Coleman-Weinberg potential at the tree-level minimum.
std::vector< double > BSMPT::Class_Potential_Origin::WeinbergThirdDerivative | ( | ) | const |
Calculates the third derivative of the Coleman-Weinberg potential at the tree-level minimum.
|
pure virtual |
This will produce a terminal output of all the model parameters. This has to be specified in the model file.
Implemented in BSMPT::Models::Class_Potential_C2HDM, BSMPT::Models::Class_Potential_CPintheDark, BSMPT::Models::Class_CxSM, BSMPT::Models::Class_Potential_N2HDM, BSMPT::Models::Class_Potential_R2HDM, BSMPT::Models::Class_SM, and BSMPT::Models::Class_Template.
|
protected |
G^{abij}
|
protected |
L_{(S),CT}^{i}
|
protected |
L_{(S),CT}^{ij}
|
protected |
L_{(S),CT}^{ijk}
|
protected |
L_{(S),CT}^{ijkl}
|
protected |
L_{(S)}^{i}
|
protected |
L_{(S)}^{ij}
|
protected |
L_{(S)}^{ijk}
|
protected |
L_{(S)}^{ijkl}
|
protected |
Y^{IJ} for Leptons
|
protected |
Y^{IJk} for Leptons
|
protected |
Y^{IJ} for Quarks
|
protected |
Y^{IJk} for Quarks
|
protected |
Storage of the Higgs rotation matrix for the Higgs mass matrix at the tree-level Vacuum
|
protected |
Storage of the Tree-Level Higgs VEV
|
protected |
Parameter to differentiate the models
|
protected |
Number of charged Higgs bosons
|
protected |
Number of gauge bosons. Do not change this in the current version as we only investigate extended Higgs sectors. If you want to extend the other sectors as well the Debye corrections have to be calculated by hand.
|
protected |
Number of all Higgs particles. This will define the size of your Higgs mass matrix.
|
protected |
Number of leptons. Do not change this in the current version as we only investigate extended Higgs sectors. If you want to extend the other sectors as well the Debye corrections have to be calculated by hand.
|
protected |
Number of neutral Higgs bosons
|
protected |
Number of Lagrange parameters in the Higgs Tree-Level potential
|
protected |
Number of Lagrange parameters in the Higgs Counterterm potential
|
protected |
Number of quarks. Do not change this in the current version as we only investigate extended Higgs sectors. If you want to extend the other sectors as well the Debye corrections have to be calculated by hand.
|
protected |
Number of VEVs you want to minimize.
|
protected |
Vector to store the parameter of the Counterterm Potential
|
protected |
Vector to store the parameter of the Potential
|
protected |
MSBar renormalization scale
std::vector<std::vector<double> > BSMPT::Class_Potential_Origin::SignSymmetries |
Find all possible sign combinations of the vevs under which the potential is invariant
|
protected |
Storage of the contributions of the counterterm potential to the triple Higgs couplings in the mass basis
|
protected |
Storage of the contributions of the Coleman-Weinberg potential to the triple Higgs couplings in the gauge basis
|
protected |
Storage of the contributions of the Coleman-Weinberg potential to the triple Higgs couplings in the mass basis
|
protected |
Storage of the contributions of the tree-level potential to the triple Higgs couplings in the mass basis
|
protected |
Variable to check if the input file has an index column or not
|
protected |
Storage of the symmetric VEV in all Higgs fields
|
protected |
Storage of the tree-level VEV in all Higgs fields
|
protected |
Storage of the tree-Level VEV in the configuration for the minimizer