/home/runner/work/kynema/kynema/kynema/src/solver/contribute_constraints_system_residual_to_vector.hpp Source File

Kynema API: /home/runner/work/kynema/kynema/kynema/src/solver/contribute_constraints_system_residual_to_vector.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
contribute_constraints_system_residual_to_vector.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Kokkos_Core.hpp>
4
5namespace kynema::solver {
6
11template <typename DeviceType>
13 template <typename ValueType>
14 using View = Kokkos::View<ValueType, DeviceType>;
15 template <typename ValueType>
17 template <typename ValueType>
18 using LeftView = Kokkos::View<ValueType, Kokkos::LayoutLeft, DeviceType>;
19
24
26 void operator()(size_t i_constraint) const {
27 const auto num_dofs = target_active_dofs(i_constraint);
28 constexpr auto force_atomic =
29 !std::is_same_v<typename DeviceType::execution_space, Kokkos::Serial>;
30 for (auto component = 0U; component < num_dofs; ++component) {
31 if constexpr (force_atomic) {
32 Kokkos::atomic_add(
35 );
36 } else {
39 }
40 }
41 }
42};
43
44} // namespace kynema::solver
Definition calculate_error_sum_squares.hpp:5
A Kernel which sums the system residual contributions for a constraint's target node into the correct...
Definition contribute_constraints_system_residual_to_vector.hpp:12
LeftView< double *[1]> residual
Definition contribute_constraints_system_residual_to_vector.hpp:23
typename View< ValueType >::const_type ConstView
Definition contribute_constraints_system_residual_to_vector.hpp:16
ConstView< double *[6]> system_residual_terms
Definition contribute_constraints_system_residual_to_vector.hpp:22
Kokkos::View< ValueType, DeviceType > View
Definition contribute_constraints_system_residual_to_vector.hpp:14
Kokkos::View< ValueType, Kokkos::LayoutLeft, DeviceType > LeftView
Definition contribute_constraints_system_residual_to_vector.hpp:18
ConstView< size_t *[6]> target_node_freedom_table
Definition contribute_constraints_system_residual_to_vector.hpp:20
KOKKOS_FUNCTION void operator()(size_t i_constraint) const
Definition contribute_constraints_system_residual_to_vector.hpp:26
ConstView< size_t * > target_active_dofs
Definition contribute_constraints_system_residual_to_vector.hpp:21