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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/math/project_points_to_target_polynomial.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
project_points_to_target_polynomial.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <array>
4#include <span>
5#include <vector>
6
7#include "interpolation.hpp"
9
10namespace kynema::math {
11
31inline std::vector<std::array<double, 3>> ProjectPointsToTargetPolynomial(
32 size_t num_inputs, size_t num_outputs, std::span<const std::array<double, 3>> input_points
33) {
34 // Calculate matrix of num_inputs points-based LSFE shape function values at ouput
35 // locations
36 const auto shape_functions = ComputeShapeFunctionValues(
37 GenerateGLLPoints(num_outputs - 1), GenerateGLLPoints(num_inputs - 1)
38 );
39
40 // Project input_points to output locations using LSFE shape functions
41 auto output_points = std::vector<std::array<double, 3>>(num_outputs);
42 for (auto output : std::views::iota(0U, num_outputs)) {
43 for (auto input : std::views::iota(0U, num_inputs)) {
44 for (auto dim : std::views::iota(0U, 3U)) {
45 output_points[output][dim] +=
46 shape_functions[input][output] * input_points[input][dim];
47 }
48 }
49 }
50 return output_points;
51}
52
53} // namespace kynema::math
Definition gll_quadrature.hpp:7
std::vector< std::array< double, 3 > > ProjectPointsToTargetPolynomial(size_t num_inputs, size_t num_outputs, std::span< const std::array< double, 3 > > input_points)
Projects 3D points from a given (lower) polynomial representation to a target (higher) polynomial rep...
Definition project_points_to_target_polynomial.hpp:31
std::vector< double > GenerateGLLPoints(const size_t order)
Generates Gauss-Lobatto-Legendre (GLL) points for spectral element discretization.
Definition interpolation.hpp:169
std::vector< std::vector< double > > ComputeShapeFunctionValues(std::span< const double > input_points, std::span< const double > output_points)
Computes shape function matrices ϕg relating points ξb to ξg At least two input points are required a...
Definition least_squares_fit.hpp:65