/home/runner/work/kynema-fmb/kynema-fmb/kynema-fmb/src/elements/beams/beam_quadrature.hpp Source File

Kynema API: /home/runner/work/kynema-fmb/kynema-fmb/kynema-fmb/src/elements/beams/beam_quadrature.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
beam_quadrature.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <algorithm>
4#include <array>
5#include <cmath>
6#include <iostream>
7#include <ranges>
8#include <span>
9#include <vector>
10
13
15
26inline std::vector<std::array<double, 2>> CreateTrapezoidalQuadrature(std::span<const double> grid) {
27 const auto n{grid.size()};
28 const auto [grid_min, grid_max] = std::ranges::minmax(grid);
29 const auto grid_range{grid_max - grid_min};
30 auto quadrature = std::vector<std::array<double, 2>>{
31 {-1., (grid[1] - grid[0]) / grid_range},
32 };
33 std::ranges::transform(
34 std::views::iota(1U, n - 1), std::back_inserter(quadrature),
35 [grid, gm = grid_min, grid_range](auto i) {
36 return std::array{
37 (2. * (grid[i] - gm) / grid_range) - 1., (grid[i + 1] - grid[i - 1]) / grid_range
38 };
39 }
40 );
41 quadrature.push_back({1., (grid[n - 1] - grid[n - 2]) / grid_range});
42 return quadrature;
43}
44
45inline std::vector<std::array<double, 2>> CreateGaussLegendreLobattoQuadrature(
46 std::span<const double> grid, std::span<const double> original_grid, size_t order
47) {
48 const auto n{grid.size()};
49 const auto [grid_min, grid_max] = std::ranges::minmax(grid);
50 const auto grid_range{grid_max - grid_min};
51 const auto sectional_weights = math::GetGllWeights(order);
52 const auto sectional_num_nodes = sectional_weights.size();
53 const auto num_sections = (n - 1) / (sectional_num_nodes - 1);
54
55 auto quadrature = std::vector<std::array<double, 2>>{};
56 std::ranges::transform(
57 grid, std::back_inserter(quadrature), [gm = grid_min, grid_range](auto grid_location) {
58 return std::array{(2. * (grid_location - gm) / grid_range) - 1., 0.};
59 }
60 );
61
62 auto section_index = 0UL;
63 for ([[maybe_unused]] auto section : std::views::iota(0U, num_sections)) {
64 const auto section_range = original_grid[section + 1] - original_grid[section];
65 const auto weight_scaling = section_range / grid_range;
66 for (auto node : std::views::iota(0U, sectional_num_nodes)) {
67 quadrature[section_index + node][1] += sectional_weights[node] * weight_scaling;
68 }
69 section_index += sectional_num_nodes - 1UL;
70 }
71 return quadrature;
72}
73
86inline std::vector<std::array<double, 2>> CreateGaussLegendreQuadrature(
87 std::span<const double> grid, std::span<const double> original_grid, size_t order
88) {
89 const auto n{grid.size()};
90 const auto [grid_min, grid_max] = std::ranges::minmax(original_grid);
91 const auto grid_range{grid_max - grid_min};
92 const auto sectional_weights = math::GetGlWeights(order);
93 const auto sectional_num_nodes = sectional_weights.size();
94 const auto num_sections = n / sectional_num_nodes;
95
96 auto quadrature = std::vector<std::array<double, 2>>{};
97 std::ranges::transform(
98 grid, std::back_inserter(quadrature), [gm = grid_min, grid_range](auto grid_location) {
99 return std::array{(2. * (grid_location - gm) / grid_range) - 1., 0.};
100 }
101 );
102
103 auto section_index = 0UL;
104 for ([[maybe_unused]] auto section : std::views::iota(0U, num_sections)) {
105 const auto section_range = original_grid[section + 1] - original_grid[section];
106 const auto weight_scaling = section_range / grid_range;
107 for (auto node : std::views::iota(0U, sectional_num_nodes)) {
108 quadrature[section_index + node][1] += sectional_weights[node] * weight_scaling;
109 }
110 section_index += sectional_num_nodes;
111 }
112 return quadrature;
113}
114
115} // namespace kynema_fmb::beams
Definition beam_quadrature.hpp:14
std::vector< std::array< double, 2 > > CreateTrapezoidalQuadrature(std::span< const double > grid)
Creates a trapezoidal quadrature rule based on a given grid.
Definition beam_quadrature.hpp:26
std::vector< std::array< double, 2 > > CreateGaussLegendreQuadrature(std::span< const double > grid, std::span< const double > original_grid, size_t order)
Creates Gauss-Legendre (GL) quadrature points and weights on [-1, 1].
Definition beam_quadrature.hpp:86
std::vector< std::array< double, 2 > > CreateGaussLegendreLobattoQuadrature(std::span< const double > grid, std::span< const double > original_grid, size_t order)
Definition beam_quadrature.hpp:45
std::vector< double > GetGllWeights(size_t order)
Definition gll_quadrature.hpp:382
std::vector< double > GetGlWeights(size_t order)
Definition gl_quadrature.hpp:185