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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/system/beams/interpolate_to_quadrature_point_for_damping.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
interpolate_to_quadrature_point_for_damping.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 double jacobian, const ConstLeftView<double*>& shape_interp,
21 const ConstLeftView<double*>& shape_deriv, const ConstView<double* [7]>& node_u,
22 const ConstView<double* [6]>& node_u_dot, const View<double[3]>& u, const View<double[4]>& r,
23 const View<double[3]>& u_prime, const View<double[4]>& r_prime, const View<double[3]>& u_dot,
24 const View<double[3]>& omega, const View<double[3]>& u_dot_prime,
25 const View<double[3]>& omega_prime
26 ) {
27 for (auto node = 0U; node < node_u.extent(0); ++node) {
28 const auto phi = shape_interp(node);
29 const auto dphiJ = shape_deriv(node) / jacobian;
30 for (auto component = 0U; component < 3U; ++component) {
31 u(component) += node_u(node, component) * phi;
32 u_prime(component) += node_u(node, component) * dphiJ;
33 u_dot(component) += node_u_dot(node, component) * phi;
34 u_dot_prime(component) += node_u_dot(node, component) * dphiJ;
35 omega(component) += node_u_dot(node, component + 3U) * phi;
36 omega_prime(component) += node_u_dot(node, component + 3U) * dphiJ;
37 }
38 for (auto component = 0U; component < 4U; ++component) {
39 r(component) += node_u(node, component + 3U) * phi;
40 r_prime(component) += node_u(node, component + 3U) * dphiJ;
41 }
42 }
43 const auto r_length = KokkosBlas::serial_nrm2(r);
44 for (auto component = 0U; component < 4U; ++component) {
46 }
47 }
48};
49} // namespace kynema::beams
Definition beam_quadrature.hpp:14
Definition interpolate_to_quadrature_point_for_damping.hpp:9
typename LeftView< ValueType >::const_type ConstLeftView
Definition interpolate_to_quadrature_point_for_damping.hpp:17
typename View< ValueType >::const_type ConstView
Definition interpolate_to_quadrature_point_for_damping.hpp:13
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_to_quadrature_point_for_damping.hpp:11
Kokkos::View< ValueType, Kokkos::LayoutLeft, DeviceType > LeftView
Definition interpolate_to_quadrature_point_for_damping.hpp:15
static KOKKOS_FUNCTION void invoke(double jacobian, const ConstLeftView< double * > &shape_interp, const ConstLeftView< double * > &shape_deriv, const ConstView< double *[7]> &node_u, const ConstView< double *[6]> &node_u_dot, const View< double[3]> &u, const View< double[4]> &r, const View< double[3]> &u_prime, const View< double[4]> &r_prime, const View< double[3]> &u_dot, const View< double[3]> &omega, const View< double[3]> &u_dot_prime, const View< double[3]> &omega_prime)
Definition interpolate_to_quadrature_point_for_damping.hpp:19