/home/runner/work/kynema/kynema/kynema/src/elements/beams/populate_element_views.hpp Source File

Kynema API: /home/runner/work/kynema/kynema/kynema/src/elements/beams/populate_element_views.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
populate_element_views.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <ranges>
4#include <span>
5#include <vector>
6
7#include "beam_element.hpp"
9#include "model/node.hpp"
10
11namespace kynema::beams {
12
14inline void PopulateNodeX0(
15 const BeamElement& elem, std::span<const Node> nodes,
16 const Kokkos::View<double* [7], Kokkos::LayoutStride, Kokkos::HostSpace>& node_x0
17) {
18 for (auto node : std::views::iota(0U, elem.node_ids.size())) {
19 for (auto component : std::views::iota(0U, 7U)) {
20 node_x0(node, component) = nodes[elem.node_ids[node]].x0[component];
21 }
22 }
23}
24
26inline void PopulateQPWeight(
27 const BeamElement& elem,
28 const Kokkos::View<double*, Kokkos::LayoutStride, Kokkos::HostSpace>& qp_weight
29) {
30 for (auto qp : std::views::iota(0U, elem.quadrature.size())) {
31 qp_weight(qp) = elem.quadrature[qp][1];
32 }
33}
34
36inline std::vector<double> MapNodePositions(const BeamElement& elem, std::span<const Node> nodes) {
37 std::vector<double> node_xi(elem.node_ids.size());
38 std::ranges::transform(elem.node_ids, std::begin(node_xi), [&](auto node_id) {
39 return 2. * nodes[node_id].s - 1.;
40 });
41 return node_xi;
42}
43
46 const BeamElement& elem, std::span<const Node> nodes,
47 const Kokkos::View<double**, Kokkos::LayoutStride, Kokkos::HostSpace>& shape_interp
48) {
49 const auto node_xi = MapNodePositions(elem, nodes);
50
51 std::vector<double> weights;
52 for (auto qp : std::views::iota(0U, elem.quadrature.size())) {
53 auto qp_xi = elem.quadrature[qp][0];
54
55 math::LagrangePolynomialInterpWeights(qp_xi, node_xi, weights);
56 for (auto node : std::views::iota(0U, node_xi.size())) {
57 shape_interp(node, qp) = weights[node];
58 }
59 }
60}
61
64 const BeamElement& elem, std::span<const Node> nodes,
65 const Kokkos::View<double**, Kokkos::LayoutStride, Kokkos::HostSpace>& shape_deriv
66) {
67 const auto node_xi = MapNodePositions(elem, nodes);
68
69 std::vector<double> weights;
70 for (auto qp : std::views::iota(0U, elem.quadrature.size())) {
71 auto qp_xi = elem.quadrature[qp][0];
72
73 math::LagrangePolynomialDerivWeights(qp_xi, node_xi, weights);
74 for (auto node : std::views::iota(0U, node_xi.size())) {
75 shape_deriv(node, qp) = weights[node];
76 }
77 }
78}
79
81inline std::vector<double> MapSectionPositions(const BeamElement& elem) {
82 std::vector<double> section_xi(elem.sections.size());
83 std::ranges::transform(elem.sections, std::begin(section_xi), [](const auto& section) {
84 return 2. * section.position - 1.;
85 });
86 return section_xi;
87}
88
90inline void PopulateQPMstar(
91 const BeamElement& elem,
92 const Kokkos::View<double* [6][6], Kokkos::LayoutStride, Kokkos::HostSpace>& qp_Mstar
93) {
94 const auto section_xi = MapSectionPositions(elem);
95 std::vector<double> weights(elem.sections.size());
96
97 for (auto qp : std::views::iota(0U, elem.quadrature.size())) {
98 auto qp_xi = elem.quadrature[qp][0];
99 math::LinearInterpWeights(qp_xi, section_xi, weights);
100 for (auto component_1 : std::views::iota(0U, 6U)) {
101 for (auto component_2 : std::views::iota(0U, 6U)) {
102 qp_Mstar(qp, component_1, component_2) = 0.;
103 }
104 }
105 for (auto section : std::views::iota(0U, section_xi.size())) {
106 for (auto component_1 : std::views::iota(0U, 6U)) {
107 for (auto component_2 : std::views::iota(0U, 6U)) {
108 qp_Mstar(qp, component_1, component_2) +=
109 elem.sections[section].M_star[component_1][component_2] * weights[section];
110 }
111 }
112 }
113 }
114}
115
117inline void PopulateQPCstar(
118 const BeamElement& elem,
119 const Kokkos::View<double* [6][6], Kokkos::LayoutStride, Kokkos::HostSpace>& qp_Cstar
120) {
121 const auto section_xi = MapSectionPositions(elem);
122 std::vector<double> weights(elem.sections.size());
123
124 for (auto qp : std::views::iota(0U, elem.quadrature.size())) {
125 auto qp_xi = elem.quadrature[qp][0];
126 math::LinearInterpWeights(qp_xi, section_xi, weights);
127 for (auto component_1 : std::views::iota(0U, 6U)) {
128 for (auto component_2 : std::views::iota(0U, 6U)) {
129 qp_Cstar(qp, component_1, component_2) = 0.;
130 }
131 }
132 for (auto section : std::views::iota(0U, section_xi.size())) {
133 for (auto component_1 : std::views::iota(0U, 6U)) {
134 for (auto component_2 : std::views::iota(0U, 6U)) {
135 qp_Cstar(qp, component_1, component_2) +=
136 elem.sections[section].C_star[component_1][component_2] * weights[section];
137 }
138 }
139 }
140 }
141}
142} // namespace kynema::beams
Definition beam_quadrature.hpp:16
void PopulateNodeX0(const BeamElement &elem, std::span< const Node > nodes, const Kokkos::View< double *[7], Kokkos::LayoutStride, Kokkos::HostSpace > &node_x0)
Populate the node initial position and orientation.
Definition populate_element_views.hpp:14
void PopulateShapeFunctionValues(const BeamElement &elem, std::span< const Node > nodes, const Kokkos::View< double **, Kokkos::LayoutStride, Kokkos::HostSpace > &shape_interp)
Populate shape function values at each quadrature point.
Definition populate_element_views.hpp:45
void PopulateQPWeight(const BeamElement &elem, const Kokkos::View< double *, Kokkos::LayoutStride, Kokkos::HostSpace > &qp_weight)
Populate the integration weights at each quadrature point.
Definition populate_element_views.hpp:26
void PopulateShapeFunctionDerivatives(const BeamElement &elem, std::span< const Node > nodes, const Kokkos::View< double **, Kokkos::LayoutStride, Kokkos::HostSpace > &shape_deriv)
Populate shape function derivatives at each quadrature point.
Definition populate_element_views.hpp:63
std::vector< double > MapSectionPositions(const BeamElement &elem)
Map section positions from [0,1] to [-1,1].
Definition populate_element_views.hpp:81
void PopulateQPCstar(const BeamElement &elem, const Kokkos::View< double *[6][6], Kokkos::LayoutStride, Kokkos::HostSpace > &qp_Cstar)
Populate stiffness matrix values at each quadrature point.
Definition populate_element_views.hpp:117
void PopulateQPMstar(const BeamElement &elem, const Kokkos::View< double *[6][6], Kokkos::LayoutStride, Kokkos::HostSpace > &qp_Mstar)
Populate mass matrix values at each quadrature point.
Definition populate_element_views.hpp:90
std::vector< double > MapNodePositions(const BeamElement &elem, std::span< const Node > nodes)
Map node positions from [0,1] to [-1,1].
Definition populate_element_views.hpp:36
void LagrangePolynomialInterpWeights(double x, std::span< const double > xs, std::vector< double > &weights)
Computes weights for Lagrange polynomial interpolation.
Definition interpolation.hpp:69
void LinearInterpWeights(double x, std::span< const double > xs, std::vector< double > &weights)
Computes weights for linear interpolation.
Definition interpolation.hpp:21
void LagrangePolynomialDerivWeights(double x, std::span< const double > xs, std::vector< double > &weights)
Computes weights for Lagrange polynomial derivative interpolation.
Definition interpolation.hpp:97
Beam element constitutes flexible beams material behavior in kynema.
Definition beam_element.hpp:17
std::vector< BeamSection > sections
Definition beam_element.hpp:20
std::vector< size_t > node_ids
Definition beam_element.hpp:19
std::vector< std::array< double, 2 > > quadrature
Definition beam_element.hpp:21