blob: dde7c708bab2e6b4bbaf0ecea8a2e4c64a53bf7f [file] [log] [blame]
Austin Schuhb0bfaf82024-06-19 19:47:23 -07001#ifndef FRC971_CONTROL_LOOPS_RUNGE_KUTTA_HELPERS_H_
2#define FRC971_CONTROL_LOOPS_RUNGE_KUTTA_HELPERS_H_
3
Austin Schuh99f7c6a2024-06-25 22:07:44 -07004#include "absl/log/check.h"
5#include "absl/log/log.h"
Austin Schuhb0bfaf82024-06-19 19:47:23 -07006#include <Eigen/Dense>
7
8namespace frc971::control_loops {
9
10// Returns a reasonable Runge Kutta initial step size. This is translated from
11// scipy.
12template <typename F, typename T>
13double SelectRungeKuttaInitialStep(const F &fn, size_t t0, T y0, T f0,
14 int error_estimator_order, double rtol,
15 double atol) {
16 constexpr int states = y0.rows();
17 const Eigen::Matrix<double, states, 1> scale =
18 atol + (y0.cwiseAbs().matrix() * rtol).array();
19 const double sqrt_rows = std::sqrt(static_cast<double>(states));
20 const double d0 = (y0.array() / scale.array()).matrix().norm() / sqrt_rows;
21 const double d1 = (f0.array() / scale.array()).matrix().norm() / sqrt_rows;
22 double h0;
23 if (d0 < 1e-5 || d1 < 1e-5) {
24 h0 = 1e-6;
25 } else {
26 h0 = 0.01 * d0 / d1;
27 }
28
29 const Eigen::Matrix<double, states, 1> y1 = y0 + h0 * f0;
30 const Eigen::Matrix<double, states, 1> f1 = fn(t0 + h0, y1);
31 const double d2 =
32 ((f1 - f0).array() / scale.array()).matrix().norm() / sqrt_rows / h0;
33
34 double h1;
35 if (d1 <= 1e-15 && d2 <= 1e-15) {
36 h1 = std::max(1e-6, h0 * 1e-3);
37 } else {
38 h1 = std::pow((0.01 / std::max(d1, d2)),
39 (1.0 / (error_estimator_order + 1.0)));
40 }
41
42 return std::min(100 * h0, h1);
43}
44
45// Performs a single step of Runge Kutta integration for the adaptive algorithm
46// below. This is translated from scipy.
47template <size_t N, size_t NStages, size_t Order, typename F>
48std::tuple<Eigen::Matrix<double, N, 1>, Eigen::Matrix<double, N, 1>> RKStep(
49 const F &fn, const double t, const Eigen::Matrix<double, N, 1> &y0,
50 const Eigen::Matrix<double, N, 1> &f0, const double h,
51 const Eigen::Matrix<double, NStages, Order> &A,
52 const Eigen::Matrix<double, 1, NStages> &B,
53 const Eigen::Matrix<double, 1, NStages> &C,
54 Eigen::Matrix<double, NStages + 1, N> &K) {
55 K.template block<N, 1>(0, 0) = f0;
56 for (size_t s = 1; s < NStages; ++s) {
57 Eigen::Matrix<double, N, 1> dy =
58 K.block(0, 0, s, N).transpose() * A.block(s, 0, 1, s).transpose() * h;
59 K.template block<1, N>(s, 0) = fn(t + C(0, s) * h, y0 + dy).transpose();
60 }
61
62 Eigen::Matrix<double, N, 1> y_new =
63 y0 + h * (K.template block<NStages, N>(0, 0).transpose() * B.transpose());
64 Eigen::Matrix<double, N, 1> f_new = fn(t + h, y_new);
65
66 K.template block<1, N>(NStages, 0) = f_new.transpose();
67
68 return std::make_tuple(y_new, f_new);
69}
70
71} // namespace frc971::control_loops
72
73#endif // FRC971_CONTROL_LOOPS_RUNGE_KUTTA_HELPERS_H_