flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
lattice_numerical_constraint.hpp
Go to the documentation of this file.
1#ifndef lattice_numerical_constraint_hpp
2#define lattice_numerical_constraint_hpp
3
4
5#include <gsl/gsl_deriv.h>
6
8
9namespace flexiblesusy {
10
12 static constexpr Real default_epsilon = 1e-8;
13};
14
17 public ForeignConstraint {
18public:
20 (std::vector<size_t> dependence, Real epsilon = default_epsilon) :
22 F_gsl.function = &c_wrap;
23 F_gsl.params = this;
24 }
25 void init(RGFlow<Lattice> *flow, size_t theory, size_t site);
26 void operator()();
28protected:
29 virtual Real c(const Real *x) const = 0;
30private:
31 std::vector<size_t> nonzeros;
32 std::vector<bool> depends_on;
33 gsl_function F_gsl;
34 static double c_wrap(double xj, void *params);
35 size_t j;
37};
38
40public:
42 (std::vector<size_t> dependence,
43 std::function<Real(const AnyNumericalConstraint *, const Real *x)> fxn,
45 NumericalConstraint(dependence, epsilon), fxn_(fxn)
46 {}
47protected:
48 Real c(const Real *x) const { return fxn_(this, x); }
49private:
50 std::function<Real(const AnyNumericalConstraint *, const Real *x)> fxn_;
51};
52
55 public ForeignMatching {
56public:
57 NumericalMatching(std::vector<size_t> depL, std::vector<size_t> depH,
59 ForeignMatching(1), nonzerosL(depL), nonzerosH(depH),
61 F_gsl.function = &c_wrap;
62 F_gsl.params = this;
63 }
64 void init(RGFlow<Lattice> *flow, size_t lower_theory);
65 void operator()();
67protected:
68 virtual Real c(const Real *w, const Real *x) const = 0;
69private:
70 std::vector<size_t> nonzerosL, nonzerosH;
71 std::vector<bool> depends_on;
72 gsl_function F_gsl;
73 static double c_wrap(double wxj, void *params);
74 size_t j;
75 Real &wx(size_t j) { return j < w.size() ? w[j] : x[j - w.size()]; }
77};
78
80public:
82 (std::vector<size_t> depL, std::vector<size_t> depH,
83 std::function<Real(const AnyNumericalMatching *,
84 const Real *w, const Real *x)> fxn,
86 NumericalMatching(depL, depH, epsilon), fxn_(fxn)
87 {}
88protected:
89 Real c(const Real *w, const Real *x) const { return fxn_(this, w, x); }
90private:
91 std::function<Real(const AnyNumericalMatching *,
92 const Real *w, const Real *x)> fxn_;
93};
94
95} // namespace flexiblesusy
96
97#endif // lattice_numerical_constraint_hpp
AnyNumericalConstraint(std::vector< size_t > dependence, std::function< Real(const AnyNumericalConstraint *, const Real *x)> fxn, Real epsilon=default_epsilon)
std::function< Real(const AnyNumericalConstraint *, const Real *x)> fxn_
std::function< Real(const AnyNumericalMatching *, const Real *w, const Real *x)> fxn_
AnyNumericalMatching(std::vector< size_t > depL, std::vector< size_t > depH, std::function< Real(const AnyNumericalMatching *, const Real *w, const Real *x)> fxn, Real epsilon=default_epsilon)
Real c(const Real *w, const Real *x) const
void init(RGFlow< Lattice > *flow, size_t theory, size_t site)
void init(RGFlow< Lattice > *flow, size_t lower_theory)
static double c_wrap(double xj, void *params)
virtual Real c(const Real *x) const =0
void init(RGFlow< Lattice > *flow, size_t theory, size_t site)
NumericalConstraint(std::vector< size_t > dependence, Real epsilon=default_epsilon)
NumericalMatching(std::vector< size_t > depL, std::vector< size_t > depH, Real epsilon=default_epsilon)
static double c_wrap(double wxj, void *params)
void init(RGFlow< Lattice > *flow, size_t lower_theory)
virtual Real c(const Real *w, const Real *x) const =0
const Real epsilon
Definition: mathdefs.hpp:18
double Real
Definition: mathdefs.hpp:12