flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
rk.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
28#ifndef RK_H
29#define RK_H
30
31#include <algorithm>
32#include <cmath>
33
34#include "logger.hpp"
35#include "error.hpp"
36#include "config.h"
37
38namespace flexiblesusy {
39
40namespace runge_kutta {
41
42namespace {
44inline double sign(double a, double b) noexcept {
45 return b >= 0 ? std::fabs(a) : -std::fabs(a);
46}
47} // anonymous namespace
48
52template <typename ArrayType, typename Derivs>
53void rungeKuttaStep(const ArrayType& y, const ArrayType& dydx, double x,
54 double h, ArrayType& yout, ArrayType& yerr, Derivs derivs)
55{
56 const double a2 = 0.2;
57 const double a3 = 0.3;
58 const double a4 = 0.6;
59 const double a5 = 1.0;
60 const double a6 = 0.875;
61 const double b21 = 0.2;
62 const double b31 = 3.0 / 40.0;
63 const double b32 = 9.0 / 40.0;
64 const double b41 = 0.3;
65 const double b42 = -0.9;
66 const double b43 = 1.2;
67 const double b51 = -11.0 / 54.0;
68 const double b52 = 2.5;
69 const double b53 = -70.0 / 27.0;
70 const double b54 = 35.0 / 27.0;
71 const double b61 = 1631.0 / 55296.0;
72 const double b62 = 175.0 / 512.0;
73 const double b63 = 575.0 / 13824.0;
74 const double b64 = 44275.0 / 110592.0;
75 const double b65 = 253.0 / 4096.0;
76 const double c1 = 37.0 / 378.0;
77 const double c3 = 250.0 / 621.0;
78 const double c4 = 125.0 / 594.0;
79 const double c6 = 512.0 / 1771.0;
80 const double dc5 = -277.00 / 14336.0;
81 const double dc1 = c1 - 2825.0 / 27648.0;
82 const double dc3 = c3 - 18575.0 / 48384.0;
83 const double dc4 = c4 - 13525.0 / 55296.0;
84 const double dc6 = c6 - 0.25;
85
86 ArrayType ytemp = b21 * h * dydx + y;
87 const ArrayType ak2 = derivs(x + a2 * h, ytemp);
88
89 // Allowing piece-wise calculating of ytemp for speed reasons
90 ytemp = y + h * (b31 * dydx + b32 * ak2);
91 const ArrayType ak3 = derivs(x + a3 * h, ytemp);
92
93 ytemp = y + h * (b41 * dydx + b42 * ak2 + b43 * ak3);
94 const ArrayType ak4 = derivs(x+a4*h,ytemp);
95
96 ytemp = y + h * (b51 * dydx + b52 * ak2 + b53 * ak3 + b54 * ak4);
97 const ArrayType ak5 = derivs(x + a5 * h, ytemp);
98
99 ytemp = y + h * (b61 * dydx + b62 * ak2 + b63 * ak3 + b64 * ak4 + b65 * ak5);
100 const ArrayType ak6 = derivs(x + a6 * h, ytemp);
101
102 yout = y + h * (c1 * dydx + c3 * ak3 + c4 * ak4 + c6 * ak6);
103 yerr = h * (dc1 * dydx + dc3 * ak3 + dc4 * ak4 + dc5 * ak5 + dc6 * ak6);
104}
105
107template <typename ArrayType, typename Derivs>
108double odeStepper(ArrayType& y, const ArrayType& dydx, double& x, double htry,
109 double eps, const ArrayType& yscal, Derivs derivs,
110 int& max_step_dir)
111{
112 const double SAFETY = 0.9;
113 const double PGROW = -0.2;
114 const double PSHRNK = -0.25;
115 const double ERRCON = 1.89e-4;
116 const int n = y.size();
117 double errmax;
118 double h = htry;
119 ArrayType yerr(n);
120 ArrayType ytemp(n);
121
122 for (;;) {
123 rungeKuttaStep(y, dydx, x, h, ytemp, yerr, derivs);
124 errmax = (yerr / yscal).abs().maxCoeff(&max_step_dir);
125 errmax /= eps;
126 if (!std::isfinite(errmax)) {
127#ifdef ENABLE_VERBOSE
128 ERROR("odeStepper: non-perturbative running at Q = "
129 << std::exp(x) << " GeV of parameter y(" << max_step_dir
130 << ") = " << y(max_step_dir) << ", dy(" << max_step_dir
131 << ")/dx = " << dydx(max_step_dir));
132#endif
133 throw NonPerturbativeRunningError(std::exp(x), max_step_dir, y(max_step_dir));
134 }
135 if (errmax <= 1.0) {
136 break;
137 }
138 const double htemp = SAFETY * h * std::pow(errmax, PSHRNK);
139 h = (h >= 0.0 ? std::max(htemp, 0.1 * h) : std::min(htemp, 0.1 * h));
140 if (x + h == x) {
141#ifdef ENABLE_VERBOSE
142 ERROR("At Q = " << std::exp(x) << " GeV "
143 "stepsize underflow in odeStepper in parameter y("
144 << max_step_dir << ") = " << y(max_step_dir) << ", dy("
145 << max_step_dir << ")/dx = " << dydx(max_step_dir));
146#endif
147 throw NonPerturbativeRunningError(std::exp(x), max_step_dir, y(max_step_dir));
148 }
149 }
150 x += h;
151 y = ytemp;
152
153 return errmax > ERRCON ? SAFETY * h * std::pow(errmax,PGROW) : 5.0 * h;
154}
155
157template <typename ArrayType, typename Derivs,
158 typename Stepper = decltype(runge_kutta::odeStepper<ArrayType,Derivs>)>
159void integrateOdes(ArrayType& ystart, double from, double to, double eps,
160 double h1, double hmin, Derivs derivs,
161 Stepper rkqs = runge_kutta::odeStepper<ArrayType,Derivs>, int max_steps = 400)
162{
163 const int nvar = ystart.size();
164 const double TINY = 1.0e-16;
165 double x = from;
166 double h = sign(h1, to - from);
167 ArrayType yscal(nvar);
168 ArrayType y(ystart);
169 ArrayType dydx;
170 int max_step_dir = 0;
171
172 for (int nstp = 0; nstp < max_steps; ++nstp) {
173 dydx = derivs(x, y);
174 yscal = y.abs() + (dydx * h).abs() + TINY;
175 if ((x + h - to) * (x + h - from) > 0.0) {
176 h = to - x;
177 }
178
179 const double hnext = rkqs(y, dydx, x, h, eps, yscal, derivs, max_step_dir);
180
181 if ((x - to) * (to - from) >= 0.0) {
182 ystart = y;
183 return;
184 }
185
186 h = hnext;
187
188 if (std::fabs(hnext) <= hmin) {
189#ifdef ENABLE_VERBOSE
190 ERROR("step size too small: |hnext| = " << hnext << " (<= " << hmin << ")");
191#endif
192 break;
193 }
194 }
195
196#ifdef ENABLE_VERBOSE
197 ERROR("Bailed out of rk.cpp:too many steps (> " << max_steps << ") in integrateOdes\n"
198 "********** Q = " << std::exp(x) << " *********");
199 ERROR("max step in direction of " << max_step_dir);
200 for (int i = 0; i < nvar; i++)
201 ERROR("y(" << i << ") = " << y(i) << " dydx(" << i <<
202 ") = " << dydx(i));
203#endif
204
205 throw NonPerturbativeRunningError(std::exp(x), max_step_dir, y(max_step_dir));
206}
207
208} // namespace runge_kutta
209
210} // namespace flexiblesusy
211
212#endif // RK_H
Non-perturbative RG running.
Definition error.hpp:142
#define ERROR(msg)
Definition logger.hpp:65
void integrateOdes(ArrayType &ystart, double from, double to, double eps, double h1, double hmin, Derivs derivs, Stepper rkqs=runge_kutta::odeStepper< ArrayType, Derivs >, int max_steps=400)
Organises integration of 1st order system of ODEs.
Definition rk.hpp:159
double odeStepper(ArrayType &y, const ArrayType &dydx, double &x, double htry, double eps, const ArrayType &yscal, Derivs derivs, int &max_step_dir)
organises the variable step-size for Runge-Kutta evolution
Definition rk.hpp:108
void rungeKuttaStep(const ArrayType &y, const ArrayType &dydx, double x, double h, ArrayType &yout, ArrayType &yerr, Derivs derivs)
Definition rk.hpp:53
T sign(T x)
Definition mathdefs.hpp:65