flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
mixings.hpp
Go to the documentation of this file.
1// ====================================================================
2// This file is part of FlexibleSUSY.
3//
4// FlexibleSUSY is free software: you can redistribute it and/or modify
5// it under the terms of the GNU General Public License as published
6// by the Free Software Foundation, either version 3 of the License,
7// or (at your option) any later version.
8//
9// FlexibleSUSY is distributed in the hope that it will be useful, but
10// WITHOUT ANY WARRANTY; without even the implied warranty of
11// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12// General Public License for more details.
13//
14// You should have received a copy of the GNU General Public License
15// along with FlexibleSUSY. If not, see
16// <http://www.gnu.org/licenses/>.
17// ====================================================================
18
19#ifndef MIXINGS_H
20#define MIXINGS_H
21
22#include "numerics2.hpp"
23#include <Eigen/Core>
24
25namespace flexiblesusy {
26
28 Eigen::Matrix<double, 1, 1>&);
29
30template<int N>
31void convert_symmetric_fermion_mixings_to_slha(Eigen::Array<double, N, 1>&,
32 Eigen::Matrix<double, N, N>&)
33{
34}
35
37 Eigen::Matrix<std::complex<double>, 1, 1>&);
38
53template<int N>
54void convert_symmetric_fermion_mixings_to_slha(Eigen::Array<double, N, 1>& m,
55 Eigen::Matrix<std::complex<double>, N, N>& z)
56{
57 for (int i = 0; i < N; i++) {
58 // check if i'th row contains non-zero imaginary parts
59 if (!is_zero(z.row(i).imag().cwiseAbs().maxCoeff())) {
60 z.row(i) *= std::complex<double>(0.0,1.0);
61 m(i) *= -1;
62 }
63 }
64}
65
67 Eigen::Matrix<double, 1, 1>&);
68
69template<int N>
70void convert_symmetric_fermion_mixings_to_hk(Eigen::Array<double, N, 1>&,
71 Eigen::Matrix<double, N, N>&)
72{
73}
74
76 Eigen::Matrix<std::complex<double>, 1, 1>&);
77
86template<int N>
87void convert_symmetric_fermion_mixings_to_hk(Eigen::Array<double, N, 1>& m,
88 Eigen::Matrix<std::complex<double>, N, N>& z)
89{
90 for (int i = 0; i < N; i++) {
91 if (m(i) < 0.) {
92 z.row(i) *= std::complex<double>(0.0,1.0);
93 m(i) *= -1;
94 }
95 }
96}
97
98} // namespace flexiblesusy
99
100#endif
std::enable_if_t< std::is_unsigned< T >::value, bool > is_zero(T a, T prec=std::numeric_limits< T >::epsilon()) noexcept
compares a number for being close to zero
Definition: numerics2.hpp:46
void convert_symmetric_fermion_mixings_to_slha(double &, Eigen::Matrix< double, 1, 1 > &)
Definition: mixings.cpp:26
void convert_symmetric_fermion_mixings_to_hk(double &, Eigen::Matrix< double, 1, 1 > &)
Definition: mixings.cpp:52