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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/elements/beams/interpolate_QP_state.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
interpolate_QP_state.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
4
6
7namespace kynema::beams {
8
12template <typename DeviceType>
14 template <typename ValueType>
15 using View = Kokkos::View<ValueType, DeviceType>;
16 template <typename ValueType>
18
19 size_t element;
20 size_t num_nodes;
24
26 void operator()(size_t qp) const {
27 auto local_total = Kokkos::Array<double, 3>{};
28 for (auto node = 0U; node < num_nodes; ++node) {
29 const auto phi = shape_interp(element, node, qp);
30 for (auto component = 0U; component < 3U; ++component) {
32 }
33 }
34 for (auto component = 0U; component < 3U; ++component) {
36 }
37 }
38};
39
43template <typename DeviceType>
45 template <typename ValueType>
46 using View = Kokkos::View<ValueType, DeviceType>;
47 template <typename ValueType>
49
50 size_t element;
51 size_t num_nodes;
56
58 void operator()(size_t qp) const {
59 const auto jacobian = qp_jacobian(element, qp);
60 auto local_total = Kokkos::Array<double, 3>{};
61 for (auto node = 0U; node < num_nodes; ++node) {
62 const auto dphi = shape_deriv(element, node, qp);
63 for (auto component = 0U; component < 3U; ++component) {
65 }
66 }
67 for (auto component = 0U; component < 3U; ++component) {
69 }
70 }
71};
72
76template <typename DeviceType>
78 template <typename ValueType>
79 using View = Kokkos::View<ValueType, DeviceType>;
80 template <typename ValueType>
82
83 size_t element;
84 size_t num_nodes;
88
90 void operator()(size_t qp) const {
91 auto local_total = Kokkos::Array<double, 4>{};
92 for (auto node = 0U; node < num_nodes; ++node) {
93 const auto phi = shape_interp(element, node, qp);
94 for (auto component = 0U; component < 4U; ++component) {
96 }
97 }
98
99 for (auto component = 0U; component < 4U; ++component) {
101 }
102 }
103};
104
108template <typename DeviceType>
110 template <typename ValueType>
111 using View = Kokkos::View<ValueType, DeviceType>;
112 template <typename ValueType>
114
115 size_t element;
116 size_t num_nodes;
121
123 void operator()(size_t qp) const {
124 const auto jacobian = qp_jacobian(element, qp);
125 auto local_total = Kokkos::Array<double, 4>{};
126 for (auto node = 0U; node < num_nodes; ++node) {
127 const auto dphi = shape_deriv(element, node, qp);
128 for (auto component = 0U; component < 4U; ++component) {
130 }
131 }
132 for (auto component = 0U; component < 4U; ++component) {
134 }
135 }
136};
137
138} // namespace kynema::beams
Definition beam_quadrature.hpp:16
KOKKOS_INLINE_FUNCTION Kokkos::Array< double, 4 > NormalizeQuaternion(const Kokkos::Array< double, 4 > &q)
Normalizes a quaternion to ensure it is a unit quaternion.
Definition quaternion_operations.hpp:309
Interpolates the rotation (r) part of the state at a quadrature point.
Definition interpolate_QP_state.hpp:77
KOKKOS_FUNCTION void operator()(size_t qp) const
Definition interpolate_QP_state.hpp:90
ConstView< double *** > shape_interp
Definition interpolate_QP_state.hpp:85
size_t element
Definition interpolate_QP_state.hpp:83
typename View< ValueType >::const_type ConstView
Definition interpolate_QP_state.hpp:81
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_QP_state.hpp:79
size_t num_nodes
Definition interpolate_QP_state.hpp:84
View< double **[4]> qp_r
Definition interpolate_QP_state.hpp:87
ConstView< double **[7]> node_u
Definition interpolate_QP_state.hpp:86
Interpolates the rotation derivative (r') part of the state at a quadrature point.
Definition interpolate_QP_state.hpp:109
ConstView< double **[7]> node_u
Definition interpolate_QP_state.hpp:119
size_t element
Definition interpolate_QP_state.hpp:115
View< double **[4]> qp_rprime
Definition interpolate_QP_state.hpp:120
size_t num_nodes
Definition interpolate_QP_state.hpp:116
KOKKOS_FUNCTION void operator()(size_t qp) const
Definition interpolate_QP_state.hpp:123
ConstView< double ** > qp_jacobian
Definition interpolate_QP_state.hpp:118
typename View< ValueType >::const_type ConstView
Definition interpolate_QP_state.hpp:113
ConstView< double *** > shape_deriv
Definition interpolate_QP_state.hpp:117
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_QP_state.hpp:111
Interpolates the displacement (u) part of the state at a quadrature point.
Definition interpolate_QP_state.hpp:13
typename View< ValueType >::const_type ConstView
Definition interpolate_QP_state.hpp:17
size_t num_nodes
Definition interpolate_QP_state.hpp:20
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_QP_state.hpp:15
View< double **[3]> qp_u
Definition interpolate_QP_state.hpp:23
ConstView< double **[7]> node_u
Definition interpolate_QP_state.hpp:22
ConstView< double *** > shape_interp
Definition interpolate_QP_state.hpp:21
size_t element
Definition interpolate_QP_state.hpp:19
KOKKOS_FUNCTION void operator()(size_t qp) const
Definition interpolate_QP_state.hpp:26
Interpolates the displacement derivative (u') part of the state at a quadrature point.
Definition interpolate_QP_state.hpp:44
ConstView< double *** > shape_deriv
Definition interpolate_QP_state.hpp:52
size_t element
Definition interpolate_QP_state.hpp:50
View< double **[3]> qp_uprime
Definition interpolate_QP_state.hpp:55
size_t num_nodes
Definition interpolate_QP_state.hpp:51
ConstView< double **[7]> node_u
Definition interpolate_QP_state.hpp:54
KOKKOS_FUNCTION void operator()(size_t qp) const
Definition interpolate_QP_state.hpp:58
typename View< ValueType >::const_type ConstView
Definition interpolate_QP_state.hpp:48
ConstView< double ** > qp_jacobian
Definition interpolate_QP_state.hpp:53
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_QP_state.hpp:46