/home/runner/work/kynema/kynema/kynema/src/math/interpolation.hpp Source File

Kynema API: /home/runner/work/kynema/kynema/kynema/src/math/interpolation.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
interpolation.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <algorithm>
4#include <cmath>
5#include <numbers>
6#include <numeric>
7#include <ranges>
8#include <span>
9#include <stdexcept>
10#include <string>
11#include <vector>
12
13namespace kynema::math {
14
22inline void LinearInterpWeights(double x, std::span<const double> xs, std::vector<double>& weights) {
23 const auto n = xs.size();
24 weights.assign(n, 0.);
25
26 const auto lower = std::ranges::lower_bound(xs, x);
27 // If x is less than the first node, first weight -> 1 and our work is done
28 if (lower == xs.begin()) {
29 weights.front() = 1.;
30 return;
31 }
32 // If x is greater than the last node, last weight -> 1 and our work is done
33 if (lower == xs.end()) {
34 weights.back() = 1.;
35 return;
36 }
37 // x is between two nodes, so compute weights for closest nodes
38 const auto index = static_cast<unsigned>(std::distance(xs.begin(), lower));
39 const auto lower_loc = xs[index - 1];
40 const auto upper_loc = xs[index];
41 const auto weight_upper = (x - lower_loc) / (upper_loc - lower_loc);
42 weights[index - 1] = 1. - weight_upper;
43 weights[index] = weight_upper;
44}
45
54inline double LinearInterp(double x, std::span<const double> xs, std::span<const double> values) {
55 std::vector<double> weights(xs.size());
56 LinearInterpWeights(x, xs, weights);
57 return std::transform_reduce(
58 std::begin(weights), std::end(weights), std::begin(values), 0., std::plus<>(),
59 std::multiplies<>()
60 );
61}
62
71 double x, std::span<const double> xs, std::vector<double>& weights
72) {
73 const auto n = xs.size();
74 weights.assign(n, 1.);
75
76 // Pre-compute (x - xs[m]) terms to avoid repeated calculations
77 std::vector<double> x_minus_xm(n);
78 for (auto m : std::views::iota(0U, n)) {
79 x_minus_xm[m] = x - xs[m];
80 }
81
82 for (auto j : std::views::iota(0U, n)) {
83 for (auto m : std::views::iota(0U, n)) {
84 if (j != m) {
85 weights[j] *= x_minus_xm[m] / (xs[j] - xs[m]);
86 }
87 }
88 }
89}
90
99 double x, std::span<const double> xs, std::vector<double>& weights
100) {
101 const auto n = xs.size();
102 weights.assign(n, 0.);
103
104 for (auto i : std::views::iota(0U, n)) {
105 auto xi = xs[i];
106 auto weight = 0.;
107 for (auto j : std::views::iota(0U, n) | std::views::filter([i](unsigned j) {
108 return j != i;
109 })) {
110 auto range = std::views::iota(0U, n) | std::views::filter([i, j](unsigned k) {
111 return k != i && k != j;
112 }) |
113 std::views::common;
114 auto prod = std::transform_reduce(
115 std::begin(range), std::end(range), 1., std::multiplies<>(),
116 [&xs, x, xi](auto k) {
117 return (x - xs[k]) / (xi - xs[k]);
118 }
119 );
120 weight += prod / (xi - xs[j]);
121 }
122 weights[i] = weight;
123 }
124}
125
137inline double LegendrePolynomial(const size_t n, const double x) {
138 // Base cases
139 if (n == 0) {
140 return 1.;
141 }
142 if (n == 1) {
143 return x;
144 }
145
146 // Compute the nth Legendre polynomial iteratively
147 double p_n_minus_2{1.}; // P_{n-2}(x)
148 double p_n_minus_1{x}; // P_{n-1}(x)
149 double p_n{0.}; // P_n(x)
150 for (auto i : std::views::iota(2U, n + 1)) {
151 const auto i_double = static_cast<double>(i);
152 p_n = ((2. * i_double - 1.) * x * p_n_minus_1 - (i_double - 1.) * p_n_minus_2) / i_double;
153 p_n_minus_2 = p_n_minus_1;
154 p_n_minus_1 = p_n;
155 }
156 return p_n;
157}
158
170inline std::vector<double> GenerateGLLPoints(const size_t order) {
171 constexpr auto max_iterations = 1000U;
172 constexpr auto convergence_tolerance = std::numeric_limits<double>::epsilon();
173
174 if (order < 1) {
175 throw std::invalid_argument("Polynomial order must be >= 1");
176 }
177
178 const size_t n_nodes = order + 1;
179 auto gll_points = std::vector<double>(n_nodes, 0.);
180 gll_points.resize(n_nodes);
181
182 // Set the endpoints fixed at [-1, 1]
183 gll_points.front() = -1.;
184 gll_points.back() = 1.;
185
186 // Find interior GLL points (1, ..., order - 1) using Newton-Raphson iteration
187 std::vector<double> legendre_poly(n_nodes, 0.);
188 for (auto i : std::views::iota(1U, order)) {
189 // Initial guess using Chebyshev-Gauss-Lobatto nodes
190 auto x_it =
191 -std::cos(static_cast<double>(i) * std::numbers::pi / static_cast<double>(order));
192
193 bool converged{false};
194 for ([[maybe_unused]] auto iteration : std::views::iota(0U, max_iterations)) {
195 const auto x_old = x_it;
196
197 // Compute Legendre polynomials up to order n
198 for (auto k : std::views::iota(0U, n_nodes)) {
199 legendre_poly[k] = LegendrePolynomial(k, x_it);
200 }
201
202 // Newton update: x_{n+1} = x_n - f(x_n)/f'(x_n)
203 const auto numerator = (x_it * legendre_poly[n_nodes - 1]) - legendre_poly[n_nodes - 2];
204 const auto denominator = static_cast<double>(n_nodes) * legendre_poly[n_nodes - 1];
205 x_it -= numerator / denominator;
206
207 // Check for convergence
208 if (std::abs(x_it - x_old) <= convergence_tolerance) {
209 converged = true;
210 break;
211 }
212 }
213
214 if (!converged) {
215 throw std::runtime_error(
216 "Newton-Raphson iteration failed to converge for GLL point index " +
217 std::to_string(i)
218 );
219 }
220
221 gll_points[i] = x_it;
222 }
223
224 return gll_points;
225}
226
227} // namespace kynema::math
Definition gl_quadrature.hpp:8
std::vector< double > GenerateGLLPoints(const size_t order)
Generates Gauss-Lobatto-Legendre (GLL) points for spectral element discretization.
Definition interpolation.hpp:170
void LagrangePolynomialInterpWeights(double x, std::span< const double > xs, std::vector< double > &weights)
Computes weights for Lagrange polynomial interpolation.
Definition interpolation.hpp:70
double LegendrePolynomial(const size_t n, const double x)
Calculates the value of Legendre polynomial of order n at point x.
Definition interpolation.hpp:137
void LinearInterpWeights(double x, std::span< const double > xs, std::vector< double > &weights)
Computes weights for linear interpolation.
Definition interpolation.hpp:22
void LagrangePolynomialDerivWeights(double x, std::span< const double > xs, std::vector< double > &weights)
Computes weights for Lagrange polynomial derivative interpolation.
Definition interpolation.hpp:98
double LinearInterp(double x, std::span< const double > xs, std::span< const double > values)
Computes linear interpolation.
Definition interpolation.hpp:54