23 const auto n = xs.size();
24 weights.assign(n, 0.);
26 const auto lower = std::ranges::lower_bound(xs, x);
28 if (lower == xs.begin()) {
33 if (lower == xs.end()) {
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;
54inline double LinearInterp(
double x, std::span<const double> xs, std::span<const double> values) {
55 std::vector<double> weights(xs.size());
57 return std::transform_reduce(
58 std::begin(weights), std::end(weights), std::begin(values), 0., std::plus<>(),
71 double x, std::span<const double> xs, std::vector<double>& weights
73 const auto n = xs.size();
74 weights.assign(n, 1.);
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];
82 for (
auto j : std::views::iota(0U, n)) {
83 for (
auto m : std::views::iota(0U, n)) {
85 weights[j] *= x_minus_xm[m] / (xs[j] - xs[m]);
99 double x, std::span<const double> xs, std::vector<double>& weights
101 const auto n = xs.size();
102 weights.assign(n, 0.);
104 for (
auto i : std::views::iota(0U, n)) {
107 for (
auto j : std::views::iota(0U, n) | std::views::filter([i](
unsigned j) {
110 auto range = std::views::iota(0U, n) | std::views::filter([i, j](
unsigned k) {
111 return k != i && k != j;
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]);
120 weight += prod / (xi - xs[j]);
147 double p_n_minus_2{1.};
148 double p_n_minus_1{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;
171 constexpr auto max_iterations = 1000U;
172 constexpr auto convergence_tolerance = std::numeric_limits<double>::epsilon();
175 throw std::invalid_argument(
"Polynomial order must be >= 1");
178 const size_t n_nodes = order + 1;
179 auto gll_points = std::vector<double>(n_nodes, 0.);
180 gll_points.resize(n_nodes);
183 gll_points.front() = -1.;
184 gll_points.back() = 1.;
187 std::vector<double> legendre_poly(n_nodes, 0.);
188 for (
auto i : std::views::iota(1U, order)) {
191 -std::cos(
static_cast<double>(i) * std::numbers::pi /
static_cast<double>(order));
193 bool converged{
false};
194 for ([[maybe_unused]]
auto iteration : std::views::iota(0U, max_iterations)) {
195 const auto x_old = x_it;
198 for (
auto k : std::views::iota(0U, n_nodes)) {
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;
208 if (std::abs(x_it - x_old) <= convergence_tolerance) {
215 throw std::runtime_error(
216 "Newton-Raphson iteration failed to converge for GLL point index " +
221 gll_points[i] = x_it;
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