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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/elements/beams/interpolate_QP_rotation.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
interpolate_QP_rotation.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
4
6
7namespace kynema::beams {
8
13template <typename DeviceType>
15 template <typename ValueType>
16 using View = Kokkos::View<ValueType, DeviceType>;
17 template <typename ValueType>
19
22 ConstView<double***> shape_interpolation; // Num Nodes x Num Quadrature points
23 ConstView<double** [7]> node_position_rotation; // Node global position vector
24 View<double** [4]> qp_rotation; // quadrature point rotation
25
27 void operator()(int element) const {
28 using Kokkos::ALL;
29 using Kokkos::make_pair;
30 using Kokkos::subview;
31
32 const auto num_nodes = num_nodes_per_element(element);
33 const auto num_qps = num_qps_per_element(element);
34 const auto shape_interp = subview(
35 shape_interpolation, element, make_pair(size_t{0U}, num_nodes),
36 make_pair(size_t{0U}, num_qps)
37 );
38 const auto node_rot = subview(
39 node_position_rotation, element, make_pair(size_t{0U}, num_nodes), make_pair(3, 7)
40 );
41 const auto qp_rot = subview(qp_rotation, element, make_pair(size_t{0U}, num_qps), ALL);
42
43 InterpQuaternion(shape_interp, node_rot, qp_rot);
44 }
45};
46
47} // namespace kynema::beams
Definition beam_quadrature.hpp:16
KOKKOS_INLINE_FUNCTION void InterpQuaternion(const shape_matrix_type &shape_matrix, const node_type &node_v, const qp_type &qp_v)
Definition interpolation_operations.hpp:44
A Kernel which interpolates a rotation quaternion on a given element from its nodes to all of it quad...
Definition interpolate_QP_rotation.hpp:14
ConstView< size_t * > num_nodes_per_element
Definition interpolate_QP_rotation.hpp:20
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_QP_rotation.hpp:16
ConstView< size_t * > num_qps_per_element
Definition interpolate_QP_rotation.hpp:21
typename View< ValueType >::const_type ConstView
Definition interpolate_QP_rotation.hpp:18
ConstView< double *** > shape_interpolation
Definition interpolate_QP_rotation.hpp:22
View< double **[4]> qp_rotation
Definition interpolate_QP_rotation.hpp:24
KOKKOS_FUNCTION void operator()(int element) const
Definition interpolate_QP_rotation.hpp:27
ConstView< double **[7]> node_position_rotation
Definition interpolate_QP_rotation.hpp:23