/home/runner/work/kynema/kynema/kynema/src/constraints/calculate_prescribed_bc_constraint.hpp Source File

Kynema API: /home/runner/work/kynema/kynema/kynema/src/constraints/calculate_prescribed_bc_constraint.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
calculate_prescribed_bc_constraint.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <KokkosBatched_Copy_Decl.hpp>
4#include <Kokkos_Core.hpp>
5
8
9namespace kynema::constraints {
10
15template <typename DeviceType>
17 template <typename ValueType>
18 using View = Kokkos::View<ValueType, DeviceType>;
19 template <typename ValueType>
21
23 const ConstView<double[3]>& X0, const ConstView<double[7]>& inputs,
24 const ConstView<double[7]>& node_u, const View<double[6]>& residual_terms,
25 const View<double[6][6]>& target_gradient_terms
26 ) {
27 using Kokkos::Array;
28 using Kokkos::make_pair;
29 using Kokkos::subview;
30 using CopyVector = KokkosBatched::SerialCopy<KokkosBatched::Trans::NoTranspose, 1>;
31 using CopyTransposeMatrix = KokkosBatched::SerialCopy<KokkosBatched::Trans::Transpose>;
32
33 const auto u1_data = Array<double, 3>{inputs(0), inputs(1), inputs(2)};
34 const auto R1_data = Array<double, 4>{inputs(3), inputs(4), inputs(5), inputs(6)};
35 const auto u2_data = Array<double, 3>{node_u(0), node_u(1), node_u(2)};
36 const auto R2_data = Array<double, 4>{node_u(3), node_u(4), node_u(5), node_u(6)};
40 auto A_data = Array<double, 9>{};
41 auto C_data = Array<double, 9>{};
43
44 const auto u1 = ConstView<double[3]>{u1_data.data()};
45 const auto R1 = ConstView<double[4]>{R1_data.data()};
46 const auto u2 = ConstView<double[3]>{u2_data.data()};
47 const auto R2 = ConstView<double[4]>{R2_data.data()};
48 const auto R1t = View<double[4]>{R1t_data.data()};
49 const auto R1_X0 = View<double[4]>{R1_X0_data.data()};
50 const auto R2_R1t = View<double[4]>{R2_R1t_data.data()};
51 const auto A = View<double[3][3]>{A_data.data()};
52 const auto C = View<double[3][3]>{C_data.data()};
53 const auto V3 = View<double[3]>{V3_data.data()};
54
55 //----------------------------------------------------------------------
56 // Residual Vector
57 //----------------------------------------------------------------------
58
59 // Phi(0:3) = u2 + X0 - u1 - R1*X0
62 for (auto component = 0; component < 3; ++component) {
63 residual_terms(component) =
65 }
66
67 // Angular residual
68 // Phi(3:6) = axial(R2*inv(RC)*inv(R1))
72 CopyVector::invoke(V3, subview(residual_terms, make_pair(3, 6)));
73
74 //----------------------------------------------------------------------
75 // Constraint Gradient Matrix
76 //----------------------------------------------------------------------
77
78 //---------------------------------
79 // Target Node
80 //---------------------------------
81
82 // B(0:3,0:3) = I
83 for (auto component = 0; component < 3; ++component) {
84 target_gradient_terms(component, component) = 1.;
85 }
86
87 // B(3:6,3:6) = AX(R1*RC*inv(R2)) = transpose(AX(R2*inv(RC)*inv(R1)))
89 CopyTransposeMatrix::invoke(
90 A, subview(target_gradient_terms, make_pair(3, 6), make_pair(3, 6))
91 );
92 }
93};
94} // namespace kynema::constraints
Definition calculate_constraint_output.hpp:8
KOKKOS_INLINE_FUNCTION void AX_Matrix(const Matrix &A, const Matrix &AX_A)
Computes AX(A) of a square matrix.
Definition matrix_operations.hpp:19
KOKKOS_INLINE_FUNCTION void QuaternionCompose(const Quaternion1 &q1, const Quaternion2 &q2, QuaternionN &qn)
Composes (i.e. multiplies) two quaternions and stores the result in a third quaternion.
Definition quaternion_operations.hpp:204
KOKKOS_INLINE_FUNCTION void RotateVectorByQuaternion(const Quaternion &q, const View1 &v, const View2 &v_rot)
Rotates provided vector by provided unit quaternion and returns the result.
Definition quaternion_operations.hpp:120
KOKKOS_INLINE_FUNCTION void QuaternionToRotationMatrix(const Quaternion &q, const RotationMatrix &R)
Converts a 4x1 quaternion to a 3x3 rotation matrix and returns the result.
Definition quaternion_operations.hpp:20
KOKKOS_INLINE_FUNCTION void QuaternionInverse(const QuaternionInput &q_in, const QuaternionOutput &q_out)
Computes the inverse of a quaternion.
Definition quaternion_operations.hpp:172
KOKKOS_INLINE_FUNCTION void AxialVectorOfMatrix(const Matrix &m, const Vector &v)
Computes the axial vector (also known as the vector representation) of a 3x3 skew-symmetric matrix.
Definition matrix_operations.hpp:47
Kernel for calculating the residual and system gradient for a Prescribed BC constraint with six degre...
Definition calculate_prescribed_bc_constraint.hpp:16
typename View< ValueType >::const_type ConstView
Definition calculate_prescribed_bc_constraint.hpp:20
Kokkos::View< ValueType, DeviceType > View
Definition calculate_prescribed_bc_constraint.hpp:18
static KOKKOS_FUNCTION void invoke(const ConstView< double[3]> &X0, const ConstView< double[7]> &inputs, const ConstView< double[7]> &node_u, const View< double[6]> &residual_terms, const View< double[6][6]> &target_gradient_terms)
Definition calculate_prescribed_bc_constraint.hpp:22