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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/system/beams/integrate_residual_vector.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
integrate_residual_vector.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
4
5namespace kynema::beams {
6
7template <typename DeviceType>
9 template <typename ValueType>
10 using View = Kokkos::View<ValueType, DeviceType>;
11 template <typename ValueType>
13 template <typename ValueType>
14 using LeftView = Kokkos::View<ValueType, Kokkos::LayoutLeft, DeviceType>;
15 template <typename ValueType>
17
18 size_t element;
19 size_t num_qps;
31
33 void operator()(size_t node) const {
34 auto local_residual = Kokkos::Array<double, 6>{};
35 for (auto qp = 0U; qp < num_qps; ++qp) {
36 const auto weight = qp_weight_(qp);
37 const auto coeff_c = weight * shape_deriv_(node, qp);
38 const auto coeff_dig = weight * qp_jacobian_(qp) * shape_interp_(node, qp);
39 for (auto component = 0U; component < 6U; ++component) {
44 }
45 }
46 for (auto component = 0U; component < 6U; ++component) {
49 }
50 }
51};
52
53} // namespace kynema::beams
Definition beam_quadrature.hpp:16
Definition integrate_residual_vector.hpp:8
ConstView< double *[6]> qp_Fc_
Definition integrate_residual_vector.hpp:25
ConstView< double *[6]> qp_Fi_
Definition integrate_residual_vector.hpp:27
typename View< ValueType >::const_type ConstView
Definition integrate_residual_vector.hpp:12
ConstView< double * > qp_weight_
Definition integrate_residual_vector.hpp:20
ConstView< double *[6]> qp_Fd_
Definition integrate_residual_vector.hpp:26
size_t num_qps
Definition integrate_residual_vector.hpp:19
ConstView< double *[6]> qp_Fe_
Definition integrate_residual_vector.hpp:28
View< double **[6]> residual_vector_terms_
Definition integrate_residual_vector.hpp:30
Kokkos::View< ValueType, Kokkos::LayoutLeft, DeviceType > LeftView
Definition integrate_residual_vector.hpp:14
Kokkos::View< ValueType, DeviceType > View
Definition integrate_residual_vector.hpp:10
ConstView< double * > qp_jacobian_
Definition integrate_residual_vector.hpp:21
KOKKOS_FUNCTION void operator()(size_t node) const
Definition integrate_residual_vector.hpp:33
ConstLeftView< double ** > shape_deriv_
Definition integrate_residual_vector.hpp:23
size_t element
Definition integrate_residual_vector.hpp:18
ConstView< double *[6]> qp_Fg_
Definition integrate_residual_vector.hpp:29
ConstView< double *[6]> node_FX_
Definition integrate_residual_vector.hpp:24
ConstLeftView< double ** > shape_interp_
Definition integrate_residual_vector.hpp:22
typename LeftView< ValueType >::const_type ConstLeftView
Definition integrate_residual_vector.hpp:16