/home/runner/work/kynema/kynema/kynema/src/system/beams/interpolate_to_quadrature_point_for_inertia.hpp Source File

Kynema API: /home/runner/work/kynema/kynema/kynema/src/system/beams/interpolate_to_quadrature_point_for_inertia.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
interpolate_to_quadrature_point_for_inertia.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <KokkosBlas.hpp>
4#include <Kokkos_Core.hpp>
5
6namespace kynema::beams {
7
8template <typename DeviceType>
10 template <typename ValueType>
11 using View = Kokkos::View<ValueType, DeviceType>;
12 template <typename ValueType>
14 template <typename ValueType>
15 using LeftView = Kokkos::View<ValueType, Kokkos::LayoutLeft, DeviceType>;
16 template <typename ValueType>
18
20 const ConstLeftView<double*>& shape_interp, const ConstView<double* [7]>& node_u,
21 const ConstView<double* [6]>& node_u_dot, const ConstView<double* [6]>& node_u_ddot,
22 const View<double[4]>& r, const View<double[3]>& u_ddot, const View<double[3]>& omega,
23 const View<double[3]>& omega_dot
24 ) {
25 for (auto node = 0U; node < node_u.extent(0); ++node) {
26 const auto phi = shape_interp(node);
27 for (auto component = 0U; component < 3U; ++component) {
28 u_ddot(component) += node_u_ddot(node, component) * phi;
29 omega(component) += node_u_dot(node, component + 3U) * phi;
30 omega_dot(component) += node_u_ddot(node, component + 3U) * phi;
31 }
32 for (auto component = 0U; component < 4U; ++component) {
33 r(component) += node_u(node, component + 3) * phi;
34 }
35 }
36 const auto r_length = KokkosBlas::serial_nrm2(r);
37 for (auto component = 0U; component < 4U; ++component) {
39 }
40 }
41};
42} // namespace kynema::beams
Definition beam_quadrature.hpp:16
Definition interpolate_to_quadrature_point_for_inertia.hpp:9
Kokkos::View< ValueType, Kokkos::LayoutLeft, DeviceType > LeftView
Definition interpolate_to_quadrature_point_for_inertia.hpp:15
typename View< ValueType >::const_type ConstView
Definition interpolate_to_quadrature_point_for_inertia.hpp:13
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_to_quadrature_point_for_inertia.hpp:11
static KOKKOS_FUNCTION void invoke(const ConstLeftView< double * > &shape_interp, const ConstView< double *[7]> &node_u, const ConstView< double *[6]> &node_u_dot, const ConstView< double *[6]> &node_u_ddot, const View< double[4]> &r, const View< double[3]> &u_ddot, const View< double[3]> &omega, const View< double[3]> &omega_dot)
Definition interpolate_to_quadrature_point_for_inertia.hpp:19
typename LeftView< ValueType >::const_type ConstLeftView
Definition interpolate_to_quadrature_point_for_inertia.hpp:17