BSMPT 3.2.1
BSMPT - Beyond the Standard Model Phase Transitions : A C++ package for the computation of the EWPT in BSM models
Loading...
Searching...
No Matches
ClassPotentialRxSM.h
Go to the documentation of this file.
1// Copyright (C) 2018 Philipp Basler and Margarete Mühlleitner
2// SPDX-FileCopyrightText: 2021 Philipp Basler, Margarete Mühlleitner and Jonas
3// Müller
4//
5// SPDX-License-Identifier: GPL-3.0-or-later
6
12#pragma once
13
14#include <string> // for string
15#include <vector> // for vector
16
18using namespace Eigen;
19
20namespace BSMPT
21{
22namespace Models
23{
24
72{
73public:
74 Class_Potential_RxSM(const ISMConstants &smConstants);
75 virtual ~Class_Potential_RxSM();
76
77 // Choice of parameters of Lagrangian from https://arxiv.org/pdf/1512.05355
78 // Eq. (11)
79 double lambdaS, lambdaHS, vS;
80
81 // Not an input parameter; lambda is fixed via the requirement of having
82 // one of the Higgs bosons as the SM one with 125.09 GeV
83 double lambda;
84
85 // Not an input parameter; set to the SM value
86 double vH;
87
88 // Not input parameters; set through the tadpole equations
89 double msq, mSsq;
90
91 double alpha;
92
93 bool UnbrokenSingletPhase;
94
95 double dmsq, dlambda, dmSsq, dlambdaS, dlambdaHS, dT1, dT2, dT3, dT4, dT5;
96
97 std::size_t pos_Gp, pos_Gm, pos_G0, pos_h, pos_H;
98 std::size_t pos_h_SM, pos_h_H;
99
100 void ReadAndSet(const std::string &linestr,
101 std::vector<double> &par) override;
102 std::vector<std::string> addLegendCT() const override;
103 std::vector<std::string> addLegendTemp() const override;
104 std::vector<std::string> addLegendTripleCouplings() const override;
105 std::vector<std::string> addLegendVEV() const override;
106
107 void set_gen(const std::vector<double> &par) override;
108 void set_CT_Pot_Par(const std::vector<double> &par) override;
109 void write() const override;
110
111 /*
112 * RxSM interaction basis:
113 * 0 1 2 3 4
114 * rho1, eta1, psi1, zeta1, zetaS
115 */
116 const std::size_t pos_rho1 = 0, pos_eta1 = 1, pos_psi1 = 2, pos_zeta1 = 3,
117 pos_zetaS = 4;
118
126 void FindMassBasisIndices(const std::vector<double> &HiggsMasses,
127 const MatrixXd &HiggsRot);
128
129 void AdjustRotationMatrix() override;
130 void TripleHiggsCouplings() override;
131 std::vector<double> calc_CT() const override;
132
133 void SetCurvatureArrays() override;
134 bool CalculateDebyeSimplified() override;
135 bool CalculateDebyeGaugeSimplified() override;
136 double VTreeSimplified(const std::vector<double> &v) const override;
137 double VCounterSimplified(const std::vector<double> &v) const override;
138 void Debugging(const std::vector<double> &input,
139 std::vector<double> &output) const override;
140};
141
142} // namespace Models
143} // namespace BSMPT
The Class_Potential_Origin class Base class for all models. This class contains all numerical calcula...
Definition ClassPotentialOrigin.h:63
The Class_Potential_RxSM class Implementation of the RxSM according to https://arxiv....
Definition ClassPotentialRxSM.h:72
std::vector< std::string > addLegendCT() const override
Definition ClassPotentialRxSM.cpp:64
std::vector< std::string > addLegendVEV() const override
Definition ClassPotentialRxSM.cpp:125
void Debugging(const std::vector< double > &input, std::vector< double > &output) const override
Definition ClassPotentialRxSM.cpp:1120
void set_gen(const std::vector< double > &par) override
Definition ClassPotentialRxSM.cpp:168
void TripleHiggsCouplings() override
Definition ClassPotentialRxSM.cpp:701
bool CalculateDebyeSimplified() override
Definition ClassPotentialRxSM.cpp:1094
std::vector< std::string > addLegendTemp() const override
Definition ClassPotentialRxSM.cpp:81
double VCounterSimplified(const std::vector< double > &v) const override
Definition ClassPotentialRxSM.cpp:1112
void SetCurvatureArrays() override
Definition ClassPotentialRxSM.cpp:795
void write() const override
Definition ClassPotentialRxSM.cpp:319
std::vector< double > calc_CT() const override
Definition ClassPotentialRxSM.cpp:388
void ReadAndSet(const std::string &linestr, std::vector< double > &par) override
Definition ClassPotentialRxSM.cpp:134
std::vector< std::string > addLegendTripleCouplings() const override
Definition ClassPotentialRxSM.cpp:92
void AdjustRotationMatrix() override
Definition ClassPotentialRxSM.cpp:629
double VTreeSimplified(const std::vector< double > &v) const override
Definition ClassPotentialRxSM.cpp:1103
void FindMassBasisIndices(const std::vector< double > &HiggsMasses, const MatrixXd &HiggsRot)
Definition ClassPotentialRxSM.cpp:500
bool CalculateDebyeGaugeSimplified() override
Definition ClassPotentialRxSM.cpp:1099
void set_CT_Pot_Par(const std::vector< double > &par) override
Definition ClassPotentialRxSM.cpp:226
This classes calculates the Bounce action of the potential with a set temperature.
Definition CalculateEtaInterface.h:24
The ISMConstants struct containing all necessary SM constants.
Definition SMparam.h:24
Definition transition_tracer.h:154