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
BSMPT::Class_Potential_Origin Class Referenceabstract

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>

Inheritance diagram for BSMPT::Class_Potential_Origin:
BSMPT::Models::Class_CxSM BSMPT::Models::Class_Potential_C2HDM BSMPT::Models::Class_Potential_CPintheDark BSMPT::Models::Class_Potential_N2HDM BSMPT::Models::Class_Potential_R2HDM BSMPT::Models::Class_SM BSMPT::Models::Class_Template

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.
 

Detailed Description

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.

Member Function Documentation

◆ addLegendCT()

virtual std::vector< std::string > BSMPT::Class_Potential_Origin::addLegendCT ( ) const
pure virtual

◆ addLegendTemp()

virtual std::vector< std::string > BSMPT::Class_Potential_Origin::addLegendTemp ( ) const
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.

◆ addLegendTripleCouplings()

virtual std::vector< std::string > BSMPT::Class_Potential_Origin::addLegendTripleCouplings ( ) const
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.

◆ addLegendVEV()

virtual std::vector< std::string > BSMPT::Class_Potential_Origin::addLegendVEV ( ) const
pure virtual

◆ boson()

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.

Parameters
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
diff0 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

◆ boson_legacy()

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

◆ calc_CT()

virtual std::vector< double > BSMPT::Class_Potential_Origin::calc_CT ( ) const
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.

◆ CalculateDebye()

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.

Parameters
forceCalculationForces the caclulation, even if the model implements a simplified version. The simplified calculation will be skipped.

◆ CalculateDebyeGauge()

void BSMPT::Class_Potential_Origin::CalculateDebyeGauge ( )

Calculates the Debye corrections to the gauge sector. By using CalculateDebyeGaugeSimplified() the runtime can be reduced.

◆ CalculateDebyeGaugeSimplified()

virtual bool BSMPT::Class_Potential_Origin::CalculateDebyeGaugeSimplified ( )
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.

◆ CalculateDebyeSimplified()

virtual bool BSMPT::Class_Potential_Origin::CalculateDebyeSimplified ( )
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.

◆ CalculatePhysicalCouplings()

void BSMPT::Class_Potential_Origin::CalculatePhysicalCouplings ( )

Calculates all triple and quartic couplings in the physical basis

◆ CalculateRatioAlpha()

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

◆ CheckImplementation()

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

◆ CheckNLOVEV()

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.

◆ CounterTerm()

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

Parameters
vthe configuration of all VEVs at which the potential should be calculated
diff0 returns the potential and i!= 0 returns the derivative of the potential w.r.t v_i
ForceExplicitCalculationCalculate the tensors directly from the tensors even if VCounterSimplified() is given

◆ CWTerm()

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.

◆ Debugging()

virtual void BSMPT::Class_Potential_Origin::Debugging ( const std::vector< double > &  input,
std::vector< double > &  output 
) const
pure virtual

◆ EWSBVEV()

double BSMPT::Class_Potential_Origin::EWSBVEV ( const std::vector< double > &  v) const
virtual

This function calculates the EW breaking VEV from all contributing field configurations.

◆ fbase()

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.

◆ fbaseFour()

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.

◆ fbaseTri()

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.

◆ FCW()

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.

◆ fermion()

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.

Parameters
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

◆ fermion_legacy()

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

◆ FirstDerivativeOfEigenvalues()

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

Parameters
M: the original matrix
MDiff: the element-wise first derivative of the matrix M with respect to the parameter you want to consider
Returns
vector containing the mass eigenvalues and then the derivatives in the same order

◆ GaugeMassesSquared()

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

Parameters
vthe configuration of all VEVs at which the eigenvalues should be evaluated
TempThe temperature at which the Debye corrected masses should be calculated
diff0 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
Returns
Vector in which the eigenvalues m^2 of the mass matrix will be stored

◆ Get_Curvature_Gauge_G2H2()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Gauge_G2H2 ( ) const
inline

Get_Curvature_Gauge_G2H2 allows const ref access to the G^{abij} Tensor.

Returns

◆ Get_Curvature_Higgs_L2()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Higgs_L2 ( ) const
inline

Get_Curvature_Higgs_L2 allows const ref access to the L_{(S)}^{ij} Tensor.

Returns

◆ Get_Curvature_Higgs_L3()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Higgs_L3 ( ) const
inline

Get_Curvature_Higgs_L3 allows const ref access to the L_{(S)}^{ijk} Tensor.

Returns

◆ Get_Curvature_Higgs_L4()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Higgs_L4 ( ) const
inline

Get_Curvature_Higgs_L4 allows const ref access to the L_{(S)}^{ijkl} Tensor.

Returns

◆ Get_Curvature_Lepton_F2()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Lepton_F2 ( ) const
inline

Get_Curvature_Lepton_F2 allows const ref access to the Y^{IJ} Tensor for Leptons.

Returns

◆ Get_Curvature_Lepton_F2H1()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Lepton_F2H1 ( ) const
inline

Get_Curvature_Lepton_F2H1 allows const ref access to the Y^{IJk} Tensor for Leptons.

Returns

◆ Get_Curvature_Quark_F2()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Quark_F2 ( ) const
inline

Get_Curvature_Quark_F2 allows const ref access to the Y^{IJ} Tensor for Quarks.

Returns

◆ Get_Curvature_Quark_F2H1()

const auto & BSMPT::Class_Potential_Origin::Get_Curvature_Quark_F2H1 ( ) const
inline

Get_Curvature_Quark_F2H1 allows const ref access to the Y^{IJk} Tensor for Quarks.

Returns

◆ get_DebyeHiggs()

const std::vector< std::vector< double > > & BSMPT::Class_Potential_Origin::get_DebyeHiggs ( ) const
inline

get_DebyeHiggs get the Debye corrections to the Higgs mass matrix

Returns

◆ get_Model()

ModelID::ModelIDs BSMPT::Class_Potential_Origin::get_Model ( ) const
inline

get_Model

Returns
ModelID of the Model

◆ get_NChargedHiggs()

std::size_t BSMPT::Class_Potential_Origin::get_NChargedHiggs ( ) const
inline

get_NChargedHiggs

Returns
NChargedHiggs

◆ get_NGauge()

std::size_t BSMPT::Class_Potential_Origin::get_NGauge ( ) const
inline

get_NGauge

Returns
NGauge

◆ get_NHiggs()

std::size_t BSMPT::Class_Potential_Origin::get_NHiggs ( ) const
inline

get_NHiggs

Returns
NHiggs

◆ get_NLepton()

std::size_t BSMPT::Class_Potential_Origin::get_NLepton ( ) const
inline

get_NLepton

Returns
NLepton

◆ get_NNeutralHiggs()

std::size_t BSMPT::Class_Potential_Origin::get_NNeutralHiggs ( ) const
inline

get_NNeutralHiggs

Returns
NNeutralHiggs

◆ get_nPar()

std::size_t BSMPT::Class_Potential_Origin::get_nPar ( ) const
inline

get_nPar

Returns
nPar

◆ get_nParCT()

std::size_t BSMPT::Class_Potential_Origin::get_nParCT ( ) const
inline

get_nParCT

Returns
nParCT

◆ get_NQuarks()

std::size_t BSMPT::Class_Potential_Origin::get_NQuarks ( ) const
inline

get_NQuarks

Returns
NQuarks

◆ get_nVEV()

std::size_t BSMPT::Class_Potential_Origin::get_nVEV ( ) const
inline

get_nVEV

Returns
nVEV

◆ get_parCTStored()

std::vector< double > BSMPT::Class_Potential_Origin::get_parCTStored ( ) const
inline

get_parCTStored

Returns
parCTStored

◆ get_parStored()

std::vector< double > BSMPT::Class_Potential_Origin::get_parStored ( ) const
inline

get_parStored

Returns
parStored

◆ get_scale()

double BSMPT::Class_Potential_Origin::get_scale ( ) const
inline

get_scale

Returns
the MSBar renormalisation scale

◆ get_TripleHiggsCorrectionsCTPhysical()

double BSMPT::Class_Potential_Origin::get_TripleHiggsCorrectionsCTPhysical ( std::size_t  i,
std::size_t  j,
std::size_t  k 
) const
inline

get_TripleHiggsCorrectionsCTPhysical

Parameters
i
j
k
Returns
TripleHiggsCorrectionsCTPhysical[i][j][k]

◆ get_TripleHiggsCorrectionsCWPhysical()

double BSMPT::Class_Potential_Origin::get_TripleHiggsCorrectionsCWPhysical ( std::size_t  i,
std::size_t  j,
std::size_t  k 
) const
inline

get_TripleHiggsCorrectionsCWPhysical

Parameters
i
j
k
Returns
TripleHiggsCorrectionsCWPhysical[i][j][k]

◆ get_TripleHiggsCorrectionsTreePhysical()

double BSMPT::Class_Potential_Origin::get_TripleHiggsCorrectionsTreePhysical ( std::size_t  i,
std::size_t  j,
std::size_t  k 
) const
inline

get_TripleHiggsCorrectionsTreePhysical

Parameters
i
j
k
Returns
TripleHiggsCorrectionsTreePhysical[i][j][k]

◆ Get_VevOrder()

const std::vector< std::size_t > & BSMPT::Class_Potential_Origin::Get_VevOrder ( ) const
inline

Get_VevOrder allows const ref access to the VevOrder vector.

Returns

◆ get_vevTreeMin() [1/2]

std::vector< double > BSMPT::Class_Potential_Origin::get_vevTreeMin ( ) const
inline

get_vevTreeMin

Returns
vevTreeMin

◆ get_vevTreeMin() [2/2]

double BSMPT::Class_Potential_Origin::get_vevTreeMin ( const std::size_t &  k) const
inline

get_vevTreeMin

Parameters
k
Returns
vevTreeMin.at(k)

◆ GetCTIdentities()

std::vector< double > BSMPT::Class_Potential_Origin::GetCTIdentities ( ) const
virtual

GetCTIdentities.

Returns
vector with the identities required to vanish for the CT calculations. Returns an empty vector if not set otherwise.

◆ getUseIndexCol()

bool BSMPT::Class_Potential_Origin::getUseIndexCol ( )

get UseIndexCol var

◆ HessianCT()

Eigen::MatrixXd BSMPT::Class_Potential_Origin::HessianCT ( const std::vector< double > &  v) const

HessianWeinberg.

Returns

◆ HiggsMassesSquared()

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

Parameters
vthe configuration of all VEVs at which the eigenvalues should be evaluated
TempThe temperature at which the Debye corrected masses should be calculated
diff0 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
Returns
Vector in which the eigenvalues m^2 of the mass matrix will be stored

◆ HiggsMassMatrix()

MatrixXd BSMPT::Class_Potential_Origin::HiggsMassMatrix ( const std::vector< double > &  v,
double  Temp = 0,
int  diff = 0 
) const

HiggsMassMatrix calculates the Higgs mass matrix.

Parameters
vthe configuration of all VEVs at which the Mass Matrix should be evaluated
TempThe temperature at which the Debye corrected masses should be calculated
diff0 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
Returns

◆ initModel() [1/2]

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

Parameters
parParameters to define the parameter point. Same input as in set_gen()
Returns
Vector with the CT

◆ initModel() [2/2]

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

◆ initVectors()

void BSMPT::Class_Potential_Origin::initVectors ( )

Initializes all vectors needed for the calculations.

◆ LeptonMasses()

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.

Parameters
vthe configuration of all VEVs at which the eigenvalues should be evaluated
Returns
Vector in which the complex eigenvalues m of the mass matrix will be stored

◆ LeptonMassesSquared()

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

Parameters
vthe configuration of all VEVs at which the eigenvalues should be evaluated
diff0 returns the masses and i!=0 returns the derivative of m^2 w.r.t v_i
Returns
Vector in which the eigenvalues m^2 of the mass matrix will be stored

◆ LeptonMassMatrix()

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 $.

Parameters
vthe configuration of all VEVs at which the matrix should be calculated
Returns
the Mass Matrix for the Leptons of the form $ M^{IJ} = Y^{IJ} + Y^{IJk} v_k $

◆ MinimizeOrderVEV()

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

Parameters
VEVminimizervector with the VEVs, used as inputs for the minimizers
Returns
vector with all VEVs, filled with zeroes for fields which do not have a VEV

◆ NablaVCT()

Eigen::VectorXd BSMPT::Class_Potential_Origin::NablaVCT ( const std::vector< double > &  v) const

NablaVCT.

Returns

◆ Prepare_Triple()

void BSMPT::Class_Potential_Origin::Prepare_Triple ( )

Sets a tensor needed to calculate the contribution of the counterterm potential to the triple Higgs couplings.

◆ QuarkMasses()

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.

Parameters
vthe configuration of all VEVs at which the eigenvalues should be evaluated
Returns
Vector in which the complex eigenvalues m of the mass matrix will be stored

◆ QuarkMassesSquared()

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.

Parameters
vthe configuration of all VEVs at which the eigenvalues should be evaluated
diff0 returns the masses and i!=0 returns the derivative of m^2 w.r.t v_i
Returns
Vector in which the eigenvalues m^2 of the mass matrix will be stored

◆ QuarkMassMatrix()

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 $.

Parameters
vthe configuration of all VEVs at which the matrix should be calculated
Returns
the Mass Matrix for the Quarks of the form $ M^{IJ} = Y^{IJ} + Y^{IJk} v_k $

◆ ReadAndSet()

virtual void BSMPT::Class_Potential_Origin::ReadAndSet ( const std::string &  linestr,
std::vector< double > &  par 
)
pure virtual

◆ resetbools()

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.

◆ resetScale()

std::vector< double > BSMPT::Class_Potential_Origin::resetScale ( const double &  newScale)

resetScale changes the MSBar scale to newScale

Parameters
newScalethe new Scale in GeV
Returns
the CT at the new scale

◆ SecondDerivativeOfEigenvaluesNonRepeated()

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!

Parameters
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
Returns
the mass eigenvalues in the vector and then the derivatives in the same order

◆ set_All()

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()

◆ set_CT_Pot_Par()

virtual void BSMPT::Class_Potential_Origin::set_CT_Pot_Par ( const std::vector< double > &  par)
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.

◆ set_gen()

virtual void BSMPT::Class_Potential_Origin::set_gen ( const std::vector< double > &  par)
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.

◆ set_InputLineNumber()

void BSMPT::Class_Potential_Origin::set_InputLineNumber ( int  InputLineNumber_in)
inline

set_InputLineNumber

Parameters
InputLineNumber_invalue to set InputLineNumber

◆ set_scale()

void BSMPT::Class_Potential_Origin::set_scale ( double  scale_new)
inline

set_scale sets the MSBar renormalisation scale to scale_new

Parameters
scale_new

◆ SetCurvatureArrays()

virtual void BSMPT::Class_Potential_Origin::SetCurvatureArrays ( )
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.

◆ SetEWVEVZero()

void BSMPT::Class_Potential_Origin::SetEWVEVZero ( std::vector< double > &  sol) const
virtual

SetEWVEVZero Set all VEV directions in sol-vector to zero that contibute to EW VEV.

Parameters
solsolution in nVEV-space
Returns
sol with EWVEV dirs set to zero

◆ setUseIndexCol()

void BSMPT::Class_Potential_Origin::setUseIndexCol ( std::string  legend)

Sets the UseIndexCol var

◆ SetUseTreeLevel()

void BSMPT::Class_Potential_Origin::SetUseTreeLevel ( bool  val)

Set the parameter UseTreeLevel to the input

◆ sym2Dim()

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

Parameters
Tensor2Dim2-dim scalar tensor
Nk1number of first indices
Nk2number of second indices

◆ sym3Dim()

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

Parameters
Tensor3Dim4-dim scalar tensor
Nk1number of first indices
Nk2number of second indices
Nk3number of third indices

◆ sym4Dim()

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

Parameters
Tensor4Dim4-dim scalar tensor
Nk1number of first indices
Nk2number of second indices
Nk3number of third indices
Nk4number of forth indices

◆ TripleHiggsCouplings()

virtual void BSMPT::Class_Potential_Origin::TripleHiggsCouplings ( )
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.

◆ V1Loop()

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.

Parameters
diff0: No derivative, > 0 derivative w.r.t to the Higgs field; -1 : derivative w.r.t to the temperature
vthe configuration of all VEVs at which the potential should be calculated
Tempthe temperature at which the potential should be evaluated
Returns
the value of the one-loop part of the effective potential

In case of diff != 0 the mass vectors will directly have the derivatives of the masses

◆ VCounterSimplified()

virtual double BSMPT::Class_Potential_Origin::VCounterSimplified ( const std::vector< double > &  v) const
pure virtual

◆ VEff()

double BSMPT::Class_Potential_Origin::VEff ( const std::vector< double > &  v,
double  Temp = 0,
int  diff = 0,
int  Order = 1 
) const
virtual

Calculates the effective potential and its derivatives.

Parameters
vvev configuration at which the potential should be evaluated
Temptemperature at which the potential should be evaluated
diffSwitch for the derivative of the potential. Default is 0 for the value of the potential
Order0 returns the tree level potential and 1 the NLO potential. Default value is the NLO potential

◆ Vl()

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.

◆ Vsb()

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.

◆ Vsf()

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.

Parameters
MassSquaredm^2
Temptemperature
nthe order to which expand the expansion
diff0 = value; >0 derivative w.r.t.

◆ VTree()

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.

Parameters
vthe configuration of all VEVs at which the potential should be calculated
diff0 returns the potential and i!= 0 returns the derivative of the potential w.r.t v_i
ForceExplicitCalculationCalculate the tensors directly from the tensors even if VTreeSimplified() is given

◆ VTreeSimplified()

virtual double BSMPT::Class_Potential_Origin::VTreeSimplified ( const std::vector< double > &  v) const
pure virtual

◆ WeinbergFirstDerivative()

std::vector< double > BSMPT::Class_Potential_Origin::WeinbergFirstDerivative ( ) const

Calculates the first derivative of the Coleman-Weinberg potential evaluated at the tree-level minimum.

◆ WeinbergForthDerivative()

std::vector< double > BSMPT::Class_Potential_Origin::WeinbergForthDerivative ( ) const

Calculates the forth derivative of the Coleman-Weinberg potential at the tree-level minimum.

◆ WeinbergSecondDerivative()

std::vector< double > BSMPT::Class_Potential_Origin::WeinbergSecondDerivative ( ) const

Calculates the second derivative of the Coleman-Weinberg potential at the tree-level minimum.

◆ WeinbergSecondDerivativeAsMatrixXd()

Eigen::MatrixXd BSMPT::Class_Potential_Origin::WeinbergSecondDerivativeAsMatrixXd ( ) const

Calculates the second derivative of the Coleman-Weinberg potential at the tree-level minimum.

◆ WeinbergThirdDerivative()

std::vector< double > BSMPT::Class_Potential_Origin::WeinbergThirdDerivative ( ) const

Calculates the third derivative of the Coleman-Weinberg potential at the tree-level minimum.

◆ write()

virtual void BSMPT::Class_Potential_Origin::write ( ) const
pure virtual

Member Data Documentation

◆ Curvature_Gauge_G2H2

std::vector<std::vector<std::vector<std::vector<double> > > > BSMPT::Class_Potential_Origin::Curvature_Gauge_G2H2
protected

G^{abij}

◆ Curvature_Higgs_CT_L1

std::vector<double> BSMPT::Class_Potential_Origin::Curvature_Higgs_CT_L1
protected

L_{(S),CT}^{i}

◆ Curvature_Higgs_CT_L2

std::vector<std::vector<double> > BSMPT::Class_Potential_Origin::Curvature_Higgs_CT_L2
protected

L_{(S),CT}^{ij}

◆ Curvature_Higgs_CT_L3

std::vector<std::vector<std::vector<double> > > BSMPT::Class_Potential_Origin::Curvature_Higgs_CT_L3
protected

L_{(S),CT}^{ijk}

◆ Curvature_Higgs_CT_L4

std::vector<std::vector<std::vector<std::vector<double> > > > BSMPT::Class_Potential_Origin::Curvature_Higgs_CT_L4
protected

L_{(S),CT}^{ijkl}

◆ Curvature_Higgs_L1

std::vector<double> BSMPT::Class_Potential_Origin::Curvature_Higgs_L1
protected

L_{(S)}^{i}

◆ Curvature_Higgs_L2

std::vector<std::vector<double> > BSMPT::Class_Potential_Origin::Curvature_Higgs_L2
protected

L_{(S)}^{ij}

◆ Curvature_Higgs_L3

std::vector<std::vector<std::vector<double> > > BSMPT::Class_Potential_Origin::Curvature_Higgs_L3
protected

L_{(S)}^{ijk}

◆ Curvature_Higgs_L4

std::vector<std::vector<std::vector<std::vector<double> > > > BSMPT::Class_Potential_Origin::Curvature_Higgs_L4
protected

L_{(S)}^{ijkl}

◆ Curvature_Lepton_F2

std::vector<std::vector<std::complex<double> > > BSMPT::Class_Potential_Origin::Curvature_Lepton_F2
protected

Y^{IJ} for Leptons

◆ Curvature_Lepton_F2H1

std::vector<std::vector<std::vector<std::complex<double> > > > BSMPT::Class_Potential_Origin::Curvature_Lepton_F2H1
protected

Y^{IJk} for Leptons

◆ Curvature_Quark_F2

std::vector<std::vector<std::complex<double> > > BSMPT::Class_Potential_Origin::Curvature_Quark_F2
protected

Y^{IJ} for Quarks

◆ Curvature_Quark_F2H1

std::vector<std::vector<std::vector<std::complex<double> > > > BSMPT::Class_Potential_Origin::Curvature_Quark_F2H1
protected

Y^{IJk} for Quarks

◆ HiggsRotationMatrix

std::vector<std::vector<double> > BSMPT::Class_Potential_Origin::HiggsRotationMatrix
protected

Storage of the Higgs rotation matrix for the Higgs mass matrix at the tree-level Vacuum

◆ HiggsVev

std::vector<double> BSMPT::Class_Potential_Origin::HiggsVev
protected

Storage of the Tree-Level Higgs VEV

◆ Model

ModelID::ModelIDs BSMPT::Class_Potential_Origin::Model = ModelID::ModelIDs::NotSet
protected

Parameter to differentiate the models

◆ NChargedHiggs

std::size_t BSMPT::Class_Potential_Origin::NChargedHiggs = 0
protected

Number of charged Higgs bosons

◆ NGauge

std::size_t BSMPT::Class_Potential_Origin::NGauge = 4
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.

◆ NHiggs

std::size_t BSMPT::Class_Potential_Origin::NHiggs = NNeutralHiggs + NChargedHiggs
protected

Number of all Higgs particles. This will define the size of your Higgs mass matrix.

◆ NLepton

std::size_t BSMPT::Class_Potential_Origin::NLepton = 9
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.

◆ NNeutralHiggs

std::size_t BSMPT::Class_Potential_Origin::NNeutralHiggs = 0
protected

Number of neutral Higgs bosons

◆ nPar

std::size_t BSMPT::Class_Potential_Origin::nPar = 0
protected

Number of Lagrange parameters in the Higgs Tree-Level potential

◆ nParCT

std::size_t BSMPT::Class_Potential_Origin::nParCT = 0
protected

Number of Lagrange parameters in the Higgs Counterterm potential

◆ NQuarks

std::size_t BSMPT::Class_Potential_Origin::NQuarks = 12
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.

◆ nVEV

std::size_t BSMPT::Class_Potential_Origin::nVEV = 0
protected

Number of VEVs you want to minimize.

◆ parCTStored

std::vector<double> BSMPT::Class_Potential_Origin::parCTStored
protected

Vector to store the parameter of the Counterterm Potential

◆ parStored

std::vector<double> BSMPT::Class_Potential_Origin::parStored
protected

Vector to store the parameter of the Potential

◆ scale

double BSMPT::Class_Potential_Origin::scale
protected

MSBar renormalization scale

◆ SignSymmetries

std::vector<std::vector<double> > BSMPT::Class_Potential_Origin::SignSymmetries

Find all possible sign combinations of the vevs under which the potential is invariant

◆ TripleHiggsCorrectionsCTPhysical

std::vector<std::vector<std::vector<double> > > BSMPT::Class_Potential_Origin::TripleHiggsCorrectionsCTPhysical
protected

Storage of the contributions of the counterterm potential to the triple Higgs couplings in the mass basis

◆ TripleHiggsCorrectionsCW

std::vector<double> BSMPT::Class_Potential_Origin::TripleHiggsCorrectionsCW
protected

Storage of the contributions of the Coleman-Weinberg potential to the triple Higgs couplings in the gauge basis

◆ TripleHiggsCorrectionsCWPhysical

std::vector<std::vector<std::vector<double> > > BSMPT::Class_Potential_Origin::TripleHiggsCorrectionsCWPhysical
protected

Storage of the contributions of the Coleman-Weinberg potential to the triple Higgs couplings in the mass basis

◆ TripleHiggsCorrectionsTreePhysical

std::vector<std::vector<std::vector<double> > > BSMPT::Class_Potential_Origin::TripleHiggsCorrectionsTreePhysical
protected

Storage of the contributions of the tree-level potential to the triple Higgs couplings in the mass basis

◆ UseIndexCol

bool BSMPT::Class_Potential_Origin::UseIndexCol = false
protected

Variable to check if the input file has an index column or not

◆ VEVSymmetric

std::vector<double> BSMPT::Class_Potential_Origin::VEVSymmetric
protected

Storage of the symmetric VEV in all Higgs fields

◆ vevTree

std::vector<double> BSMPT::Class_Potential_Origin::vevTree
protected

Storage of the tree-level VEV in all Higgs fields

◆ vevTreeMin

std::vector<double> BSMPT::Class_Potential_Origin::vevTreeMin
protected

Storage of the tree-Level VEV in the configuration for the minimizer


The documentation for this class was generated from the following files: