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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/elements/beams/beams.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
beams.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
4
6
7namespace kynema {
8
21template <typename DeviceType>
22struct Beams {
23 template <typename ValueType>
24 using View = Kokkos::View<ValueType, DeviceType>;
25
26 size_t num_elems; // Total number of element
27 size_t max_elem_nodes; // Maximum number of nodes per element
28 size_t max_elem_qps; // Maximum number of quadrature points per element
29
30 // Element-based data
33 View<size_t**> node_state_indices; // State row index for each node
36 View<double* [6]> element_mu; // Damping coefficients per element
37
39
40 // Node-based data
41 View<double** [7]> node_x0; // Initial position/rotation
42 View<double** [7]> node_u; // State: translation/rotation displacement
43 View<double** [6]> node_u_dot; // State: translation/rotation velocity
44 View<double** [6]> node_u_ddot; // State: translation/rotation acceleration
45 View<double** [6]> node_FX; // External forces
46
47 // Quadrature point data
48 View<double**> qp_weight; // Integration weights
49 View<double**> qp_jacobian; // Jacobian vector
50 View<double** [6][6]> qp_Mstar; // Mass matrix in material frame
51 View<double** [6][6]> qp_Cstar; // Stiffness matrix in material frame
52 View<double** [7]> qp_x; // Current position/orientation
53 View<double** [3]> qp_x0; // Initial position
54 View<double** [3]> qp_x0_prime; // Initial position derivative
55 View<double** [4]> qp_r0; // Initial rotation
56 View<double** [3]> qp_u; // State: translation displacement
57 View<double** [3]> qp_u_prime; // State: translation displacement derivative
58 View<double** [3]> qp_u_dot; // State: translation velocity
59 View<double** [3]> qp_u_ddot; // State: translation acceleration
60 View<double** [4]> qp_r; // State: rotation
61 View<double** [4]> qp_r_prime; // State: rotation derivative
62 View<double** [3]> qp_omega; // State: angular velocity
63 View<double** [3]> qp_omega_dot; // State: position/rotation
64 View<double** [3]> qp_deformation; // Deformation relative to rigid body motion
65 View<double** [3][4]> qp_E; // Quaternion derivative
66 View<double** [6]> qp_Fe; // External force
67
70
71 // Shape Function data
72 View<double***> shape_interp; // Shape function values
73 View<double***> shape_deriv; // Shape function derivatives
74
75 // Constructor which initializes views based on given sizes
76 Beams(const size_t n_beams, const size_t max_e_nodes, const size_t max_e_qps)
80 // Element Data
82 Kokkos::view_alloc("num_nodes_per_element", Kokkos::WithoutInitializing), num_elems
83 ),
85 Kokkos::view_alloc("num_qps_per_element", Kokkos::WithoutInitializing), num_elems
86 ),
88 Kokkos::view_alloc("node_state_indices", Kokkos::WithoutInitializing), num_elems,
90 ),
92 Kokkos::view_alloc("element_freedom_signature", Kokkos::WithoutInitializing),
94 ),
96 Kokkos::view_alloc("element_freedom_table", Kokkos::WithoutInitializing), num_elems,
98 ),
101 // Node Data
102 node_x0(
104 ),
105 node_u(
107 ),
111 ),
115 ),
116 node_FX(
117 Kokkos::view_alloc("node_force_external", Kokkos::WithoutInitializing), num_elems,
119 ),
120 // Quadrature Point data
121 qp_weight(
123 ),
126 ),
127 qp_Mstar(
129 ),
130 qp_Cstar(
132 ),
137 ),
142 ),
143 qp_u_dot(
145 ),
146 qp_u_ddot(
148 ),
152 ),
153 qp_omega(
155 ),
159 ),
163 ),
166 Kokkos::view_alloc("residual_vector_terms", Kokkos::WithoutInitializing), num_elems,
168 ),
170 Kokkos::view_alloc("system_matrix_terms", Kokkos::WithoutInitializing), num_elems,
172 ),
173 // Shape Function data
177 ),
181 ) {
183 }
184};
185
186} // namespace kynema
Definition calculate_constraint_output.hpp:8
Contains the field variables needed to compute the per-element contributions to the residual vector a...
Definition beams.hpp:22
View< double **[3]> qp_u_ddot
Definition beams.hpp:59
View< double ** > qp_weight
Definition beams.hpp:48
size_t num_elems
Definition beams.hpp:26
View< double *** > shape_interp
Definition beams.hpp:72
View< double **[3]> qp_omega_dot
Definition beams.hpp:63
View< double **[3]> qp_deformation
Definition beams.hpp:64
View< double **[4]> qp_r
Definition beams.hpp:60
size_t max_elem_nodes
Definition beams.hpp:27
View< double **[4]> qp_r_prime
Definition beams.hpp:61
View< dof::FreedomSignature ** > element_freedom_signature
Definition beams.hpp:34
View< double **[7]> qp_x
Definition beams.hpp:52
View< double *** > shape_deriv
Definition beams.hpp:73
View< double **[6]> node_u_ddot
Definition beams.hpp:44
View< double **[3]> qp_x0_prime
Definition beams.hpp:54
View< double **[6][6]> qp_Cstar
Definition beams.hpp:51
View< size_t * > num_qps_per_element
Definition beams.hpp:32
Kokkos::View< ValueType, DeviceType > View
Definition beams.hpp:24
View< double **[4]> qp_r0
Definition beams.hpp:55
View< double **[7]> node_x0
Definition beams.hpp:41
View< double ** > qp_jacobian
Definition beams.hpp:49
View< size_t ** > node_state_indices
Definition beams.hpp:33
View< double **[6]> node_FX
Definition beams.hpp:45
View< double **[6][6]> qp_Mstar
Definition beams.hpp:50
View< size_t **[6]> element_freedom_table
Definition beams.hpp:35
View< double **[3]> qp_u
Definition beams.hpp:56
View< double[3]> gravity
Definition beams.hpp:38
View< double **[3]> qp_u_dot
Definition beams.hpp:58
Beams(const size_t n_beams, const size_t max_e_nodes, const size_t max_e_qps)
Definition beams.hpp:76
View< double **[3]> qp_omega
Definition beams.hpp:62
View< double **[6]> residual_vector_terms
Definition beams.hpp:68
View< double **[3]> qp_x0
Definition beams.hpp:53
View< double *[6]> element_mu
Definition beams.hpp:36
View< double **[6]> qp_Fe
Definition beams.hpp:66
View< double ***[6][6]> system_matrix_terms
Definition beams.hpp:69
View< size_t * > num_nodes_per_element
Definition beams.hpp:31
View< double **[6]> node_u_dot
Definition beams.hpp:43
size_t max_elem_qps
Definition beams.hpp:28
View< double **[3]> qp_u_prime
Definition beams.hpp:57
View< double **[3][4]> qp_E
Definition beams.hpp:65
View< double **[7]> node_u
Definition beams.hpp:42