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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/elements/beams/interpolate_to_quadrature_points.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
interpolate_to_quadrature_points.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
4
8
9namespace kynema::beams {
10
22template <typename DeviceType>
24 using TeamPolicy = Kokkos::TeamPolicy<typename DeviceType::execution_space>;
25 using member_type = typename TeamPolicy::member_type;
26 template <typename ValueType>
27 using View = Kokkos::View<ValueType, DeviceType>;
28 template <typename ValueType>
30
32 ConstView<size_t*> num_qps_per_element; //< Number of QPs per element
33 ConstView<double***> shape_interp; //< shape functions at QPs
34 ConstView<double***> shape_deriv; //< shape function derivatives at QPs
35 ConstView<double**> qp_jacobian; //< Jacobian at QPs
36 ConstView<double** [7]> node_u; //< Nodal displacements
37 ConstView<double** [6]> node_u_dot; //< Nodal velocities
38 ConstView<double** [6]> node_u_ddot; //< Nodal accelerations
39 ConstView<double** [3]> qp_x0; //< Initial positions at QPs
40 ConstView<double** [4]> qp_r0; //< Initial rotations at QPs
41
42 // Output quantities at quadrature points
43 View<double** [3]> qp_u; //< Interpolated displacements at QPs
44 View<double** [3]> qp_uprime; //< Displacement derivatives at QPs
45 View<double** [4]> qp_r; //< Interpolated rotations at QPs
46 View<double** [4]> qp_rprime; //< Rotation derivatives at QPs
47 View<double** [3]> qp_u_dot; //< Interpolated velocities at QPs
48 View<double** [3]> qp_omega; //< Interpolated angular velocities at QPs
49 View<double** [3]> qp_u_ddot; //< Interpolated accelerations at QPs
50 View<double** [3]> qp_omega_dot; //< Interpolated angular accelerations at QPs
51 View<double** [7]> qp_x; //< Final positions of quadrature points
52
55 using Kokkos::ALL;
56 using Kokkos::make_pair;
57 using Kokkos::parallel_for;
58
59 const auto element = static_cast<size_t>(member.league_rank());
60 const auto num_nodes = num_nodes_per_element(element);
61 const auto num_qps = num_qps_per_element(element);
62
63 auto qp_range = Kokkos::TeamThreadRange(member, num_qps);
64
68 );
72 element, num_nodes, shape_deriv, qp_jacobian, node_u, qp_uprime
73 }
74 );
78 );
82 element, num_nodes, shape_deriv, qp_jacobian, node_u, qp_rprime
83 }
84 );
88 element, num_nodes, shape_interp, subview(node_u_dot, ALL, ALL, make_pair(0, 3)),
90 }
91 );
95 element, num_nodes, shape_interp, subview(node_u_dot, ALL, ALL, make_pair(3, 6)),
97 }
98 );
100 qp_range,
102 element, num_nodes, shape_interp, subview(node_u_ddot, ALL, ALL, make_pair(0, 3)),
104 }
105 );
107 qp_range,
109 element, num_nodes, shape_interp, subview(node_u_ddot, ALL, ALL, make_pair(3, 6)),
111 }
112 );
115 );
116 }
117};
118
119} // namespace kynema::beams
Definition beam_quadrature.hpp:16
Interpolates various quantities from nodes to quadrature points for beam elements.
Definition interpolate_to_quadrature_points.hpp:23
View< double **[3]> qp_omega_dot
Definition interpolate_to_quadrature_points.hpp:50
View< double **[3]> qp_u_ddot
Definition interpolate_to_quadrature_points.hpp:49
ConstView< double ** > qp_jacobian
Definition interpolate_to_quadrature_points.hpp:35
ConstView< double **[6]> node_u_ddot
Definition interpolate_to_quadrature_points.hpp:38
typename TeamPolicy::member_type member_type
Definition interpolate_to_quadrature_points.hpp:25
ConstView< double *** > shape_interp
Definition interpolate_to_quadrature_points.hpp:33
View< double **[7]> qp_x
Definition interpolate_to_quadrature_points.hpp:51
ConstView< size_t * > num_nodes_per_element
Definition interpolate_to_quadrature_points.hpp:31
Kokkos::View< ValueType, DeviceType > View
Definition interpolate_to_quadrature_points.hpp:27
View< double **[3]> qp_u_dot
Definition interpolate_to_quadrature_points.hpp:47
View< double **[3]> qp_omega
Definition interpolate_to_quadrature_points.hpp:48
View< double **[4]> qp_r
Definition interpolate_to_quadrature_points.hpp:45
ConstView< double **[4]> qp_r0
Definition interpolate_to_quadrature_points.hpp:40
ConstView< size_t * > num_qps_per_element
Definition interpolate_to_quadrature_points.hpp:32
ConstView< double **[6]> node_u_dot
Definition interpolate_to_quadrature_points.hpp:37
ConstView< double **[3]> qp_x0
Definition interpolate_to_quadrature_points.hpp:39
View< double **[4]> qp_rprime
Definition interpolate_to_quadrature_points.hpp:46
KOKKOS_FUNCTION void operator()(member_type member) const
Definition interpolate_to_quadrature_points.hpp:54
ConstView< double *** > shape_deriv
Definition interpolate_to_quadrature_points.hpp:34
typename View< ValueType >::const_type ConstView
Definition interpolate_to_quadrature_points.hpp:29
Kokkos::TeamPolicy< typename DeviceType::execution_space > TeamPolicy
Definition interpolate_to_quadrature_points.hpp:24
View< double **[3]> qp_uprime
Definition interpolate_to_quadrature_points.hpp:44
ConstView< double **[7]> node_u
Definition interpolate_to_quadrature_points.hpp:36
View< double **[3]> qp_u
Definition interpolate_to_quadrature_points.hpp:43