Skip to content

Commit

Permalink
Add unit test for CMSSM coefficients calculation
Browse files Browse the repository at this point in the history
  • Loading branch information
Dylan Harries committed Jan 25, 2017
1 parent 049c81e commit ea88aed
Show file tree
Hide file tree
Showing 4 changed files with 224 additions and 1 deletion.
22 changes: 21 additions & 1 deletion test/module.mk
Expand Up @@ -51,7 +51,7 @@ TEST_SH := \

ifneq ($(findstring lattice,$(SOLVERS)),)
TEST_SRC +=
endif
endif # ifneq($(findstring lattice,$(SOLVERS)),)

ifneq ($(findstring two_scale,$(SOLVERS)),)
TEST_SRC += \
Expand Down Expand Up @@ -149,7 +149,25 @@ ifeq ($(WITH_sm) $(WITH_smcw),yes yes)
TEST_SRC += \
$(DIR)/test_two_scale_sm_smcw_integration.cpp
endif

endif # ifneq ($(findstring two_scale,$(SOLVERS)),)

ifneq ($(findstring semi_analytic,$(SOLVERS)),)

ifeq ($(WITH_CMSSMSemiAnalytic), yes)
TEST_SRC += \
$(DIR)/test_CMSSMSemiAnalytic_semi_analytic_solutions.cpp
endif

ifneq ($(findstring two_scale,$(SOLVERS)),)
ifeq ($(WITH_CMSSM) $(WITH_CMSSMSemiAnalytic), yes yes)
TEST_SH += \
$(DIR)/test_CMSSMSemiAnalytic_spectrum.sh
endif
endif

endif # ifneq ($(findstring semi_analytic,$(SOLVERS)),)

ifeq ($(WITH_CMSSM) $(WITH_NMSSM),yes yes)
TEST_SRC += \
$(DIR)/test_CMSSM_NMSSM_linking.cpp
Expand Down Expand Up @@ -758,6 +776,8 @@ $(DIR)/test_NSM_low_scale_constraint.x: $(LIBNSM) $(LIBFLEXI) $(filter-out -%,$(

$(DIR)/test_VCMSSM_ewsb.x: $(LIBVCMSSM) $(LIBCMSSM) $(LIBFLEXI) $(LIBLEGACY) $(filter-out -%,$(LOOPFUNCLIBS))

$(DIR)/test_CMSSMSemiAnalytic_semi_analytic_solutions.x: $(LIBCMSSMSemiAnalytic) $(LIBFLEXI) $(LIBLEGACY) $(filter-out -%,$(LOOPFUNCLIBS))

# general test rule which links all libraries needed for a generated model
$(DIR)/test_%.x: $(DIR)/test_%.o
$(CXX) -o $@ $(call abspathx,$^) $(filter -%,$(LOOPFUNCLIBS)) $(BOOSTTESTLIBS) $(BOOSTTHREADLIBS) $(THREADLIBS) $(GSLLIBS) $(LAPACKLIBS) $(BLASLIBS) $(FLIBS)
Expand Down
18 changes: 18 additions & 0 deletions test/test.h
Expand Up @@ -157,6 +157,24 @@ void check_relative_dev(T a, T b, const std::string& testMsg, T max_dev)
}
}

template <class Derived>
void check_relative_dev(const Eigen::MatrixBase<Derived>& a,
const Eigen::MatrixBase<Derived>& b,
const std::string& testMsg,
typename Derived::Scalar max_dev)
{
assert(a.rows() == b.rows());
assert(a.cols() == b.cols());

for (int i = 0; i < a.rows(); ++i) {
for (int l = 0; l < a.cols(); ++l) {
std::ostringstream element;
element << testMsg << " [element " << i << "," << l << "]";
check_relative_dev(a(i,l), b(i,l), element.str(), max_dev);
}
}
}

} // namespace flexiblesusy

#define S(x) #x
Expand Down
88 changes: 88 additions & 0 deletions test/test_CMSSMSemiAnalytic.hpp
@@ -0,0 +1,88 @@

#ifndef TEST_CMSSMSEMIANALYTIC_H
#define TEST_CMSSMSEMIANALYTIC_H

#include "CMSSMSemiAnalytic_input_parameters.hpp"
#include "CMSSMSemiAnalytic_mass_eigenstates.hpp"
#include "wrappers.hpp"
#include "ew_input.hpp"

void setup_CMSSMSemiAnalytic_const(flexiblesusy::CMSSMSemiAnalytic_mass_eigenstates& m,
const flexiblesusy::CMSSMSemiAnalytic_input_parameters& input)
{
using namespace flexiblesusy;

const double ALPHASMZ = 0.1176;
const double ALPHAMZ = 1.0 / 127.918;
const double sinthWsq = 0.23122;
const double alpha1 = 5.0 * ALPHAMZ / (3.0 * (1.0 - sinthWsq));
const double alpha2 = ALPHAMZ / sinthWsq;
const double g1 = Sqrt(4.0 * Pi * alpha1);
const double g2 = Sqrt(4.0 * Pi * alpha2);
const double g3 = Sqrt(4.0 * Pi * ALPHASMZ);
const double tanBeta = input.TanBeta;
const double sinBeta = sin(atan(tanBeta));
const double cosBeta = cos(atan(tanBeta));
const double M12 = input.m12;
const double m0 = input.m12;
const double a0 = input.Azero + 100.;
const double root2 = Sqrt(2.0);
const double vev = 246.0;
const double vu = vev * sinBeta;
const double vd = vev * cosBeta;
const double susyMu = input.MuInput;
const double BMu = Sqr(2.0 * susyMu);
const double scale = Electroweak_constants::MZ;

Eigen::Matrix<double,3,3> Yu(Eigen::Matrix<double,3,3>::Zero());
Eigen::Matrix<double,3,3> Yd(Eigen::Matrix<double,3,3>::Zero());
Eigen::Matrix<double,3,3> Ye(Eigen::Matrix<double,3,3>::Zero());
Eigen::Matrix<double,3,3> mm0(Eigen::Matrix<double,3,3>::Zero());

Yu(2,2) = 165.0 * root2 / (vev * sinBeta);
Yd(2,2) = 2.9 * root2 / (vev * cosBeta);
Ye(2,2) = 1.77699 * root2 / (vev * cosBeta);

mm0 = Sqr(m0) * Eigen::Matrix<double,3,3>::Identity();

m.set_input_parameters(input);
m.set_scale(scale);
m.set_loops(1);
m.set_thresholds(3);
m.set_g1(g1);
m.set_g2(g2);
m.set_g3(g3);
m.set_Yu(Yu);
m.set_Yd(Yd);
m.set_Ye(Ye);
m.set_MassB(M12);
m.set_MassG(M12);
m.set_MassWB(M12);
m.set_mq2(mm0);
m.set_ml2(mm0);
m.set_md2(mm0);
m.set_mu2(mm0);
m.set_me2(mm0);
m.set_mHd2(Sqr(m0));
m.set_mHu2(Sqr(m0));
m.set_TYu(a0 * Yu);
m.set_TYd(a0 * Yd);
m.set_TYe(a0 * Ye);
m.set_Mu(susyMu);
m.set_BMu(BMu);
m.set_vu(vu);
m.set_vd(vd);
}

void setup_CMSSMSemiAnalytic(flexiblesusy::CMSSMSemiAnalytic_mass_eigenstates& m,
flexiblesusy::CMSSMSemiAnalytic_input_parameters& input)
{
input.m12 = 500.;
input.Azero = 100.;
input.TanBeta = 7.0;
input.MuInput = 200.;

setup_CMSSMSemiAnalytic_const(m, input);
}

#endif
97 changes: 97 additions & 0 deletions test/test_CMSSMSemiAnalytic_semi_analytic_solutions.cpp
@@ -0,0 +1,97 @@
#define BOOST_TEST_DYN_LINK
#define BOOST_TEST_MODULE test_CMSSMSemiAnalytic_semi_analytic_solutions

#include <boost/test/unit_test.hpp>

#include "test.h"
#include "test_CMSSMSemiAnalytic.hpp"
#include "CMSSMSemiAnalytic_semi_analytic_solutions.hpp"

using namespace flexiblesusy;

struct Boundary_values {
double m12{};
double Azero{};
double m0Sq{};
double BMu0{};
double Mu{};
};

void setup_high_scale_CMSSMSemiAnalytic(
CMSSMSemiAnalytic_mass_eigenstates& model, const Boundary_values& values)
{
model.set_TYu(model.get_Yu() * values.Azero);
model.set_TYd(model.get_Yd() * values.Azero);
model.set_TYe(model.get_Ye() * values.Azero);

model.set_mq2(values.m0Sq * UNITMATRIX(3));
model.set_mu2(values.m0Sq * UNITMATRIX(3));
model.set_md2(values.m0Sq * UNITMATRIX(3));
model.set_ml2(values.m0Sq * UNITMATRIX(3));
model.set_me2(values.m0Sq * UNITMATRIX(3));
model.set_mHd2(values.m0Sq);
model.set_mHu2(values.m0Sq);

model.set_MassB(values.m12);
model.set_MassWB(values.m12);
model.set_MassG(values.m12);

model.set_Mu(values.Mu);
model.set_BMu(values.BMu0);
}

BOOST_AUTO_TEST_CASE( test_CMSSMSemiAnalytic_coefficients )
{
CMSSMSemiAnalytic_input_parameters input;
CMSSMSemiAnalytic_mass_eigenstates model(input);
setup_CMSSMSemiAnalytic(model, input);

const double high_scale = 2.e16;
model.run_to(high_scale);

Boundary_values values;
values.m12 = 300.0;
values.Azero = -200.;
values.m0Sq = Sqr(250);
values.BMu0 = -1.e4;
values.Mu = model.get_Mu();

setup_high_scale_CMSSMSemiAnalytic(model, values);

CMSSMSemiAnalytic_semi_analytic_solutions solns;
solns.set_input_scale(high_scale);
solns.set_output_scale(Electroweak_constants::MZ);
solns.set_AzeroBasis(values.Azero);
solns.set_m12Basis(values.m12);
solns.set_m0SqBasis(values.m0Sq);
solns.set_BMu0Basis(values.BMu0);
solns.set_MuBasis(values.Mu);

solns.calculate_coefficients(model);

model.run_to(Electroweak_constants::MZ);

CMSSMSemiAnalytic_mass_eigenstates coeffs_model(model);
solns.evaluate_solutions(coeffs_model);

BOOST_CHECK_CLOSE_FRACTION(model.get_MassB(), coeffs_model.get_MassB(), 1.0e-3);
BOOST_CHECK_CLOSE_FRACTION(model.get_MassWB(), coeffs_model.get_MassWB(), 1.0e-3);
BOOST_CHECK_CLOSE_FRACTION(model.get_MassG(), coeffs_model.get_MassG(), 1.0e-3);

BOOST_CHECK_CLOSE_FRACTION(model.get_BMu(), coeffs_model.get_BMu(), 1.0e-3);

BOOST_CHECK_CLOSE_FRACTION(model.get_mHd2(), coeffs_model.get_mHd2(), 1.0e-3);
BOOST_CHECK_CLOSE_FRACTION(model.get_mHu2(), coeffs_model.get_mHu2(), 1.0e-3);

TEST_CLOSE_REL(model.get_TYu(), coeffs_model.get_TYu(), 1.0e-3);
TEST_CLOSE_REL(model.get_TYd(), coeffs_model.get_TYd(), 1.0e-3);
TEST_CLOSE_REL(model.get_TYe(), coeffs_model.get_TYe(), 1.0e-3);

TEST_CLOSE_REL(model.get_mq2(), coeffs_model.get_mq2(), 1.0e-2);
TEST_CLOSE_REL(model.get_mu2(), coeffs_model.get_mu2(), 1.0e-2);
TEST_CLOSE_REL(model.get_md2(), coeffs_model.get_md2(), 1.0e-2);
TEST_CLOSE_REL(model.get_ml2(), coeffs_model.get_ml2(), 1.0e-2);
TEST_CLOSE_REL(model.get_me2(), coeffs_model.get_me2(), 1.0e-2);

BOOST_CHECK_EQUAL(get_errors(), 0);
}

0 comments on commit ea88aed

Please sign in to comment.