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

Kynema API: /home/runner/work/kynema/kynema/kynema/src/solver/contribute_masses_to_sparse_matrix.hpp Source File
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
contribute_masses_to_sparse_matrix.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <KokkosSparse.hpp>
4#include <Kokkos_Core.hpp>
5
7
8namespace kynema::solver {
9
14template <typename CrsMatrixType>
16 using DeviceType = typename CrsMatrixType::device_type;
17 using RowDataType = typename CrsMatrixType::values_type::non_const_type;
18 using ColIdxType = typename CrsMatrixType::staticcrsgraph_type::entries_type::non_const_type;
19 using TeamPolicy = Kokkos::TeamPolicy<typename DeviceType::execution_space>;
20 using member_type = typename TeamPolicy::member_type;
21 template <typename ValueType>
22 using View = Kokkos::View<ValueType, DeviceType>;
23 template <typename ValueType>
25
26 double conditioner{};
30 CrsMatrixType sparse;
31
34 const auto element = member.league_rank();
35 constexpr auto is_sorted = true;
36 constexpr auto force_atomic =
37 !std::is_same_v<typename DeviceType::execution_space, Kokkos::Serial>;
38
39 constexpr auto num_dofs = 6;
40 constexpr auto hint = static_cast<typename CrsMatrixType::ordinal_type>(0);
41
42 const auto first_column =
43 static_cast<typename CrsMatrixType::ordinal_type>(element_freedom_table(element, 0));
44 Kokkos::parallel_for(Kokkos::TeamVectorRange(member, num_dofs), [&](int component_1) {
45 const auto row_num = static_cast<int>(element_freedom_table(element, component_1));
46 auto row = sparse.row(row_num);
47 auto offset = KokkosSparse::findRelOffset(
48 &(row.colidx(0)), row.length, first_column, hint, is_sorted
49 );
50 for (auto component_2 = 0; component_2 < num_dofs; ++component_2, ++offset) {
51 const auto contribution = dense(element, component_1, component_2) * conditioner;
52 if constexpr (force_atomic) {
53 Kokkos::atomic_add(&(row.value(offset)), contribution);
54 } else {
55 row.value(offset) += contribution;
56 }
57 }
58 });
59 }
60};
61
62} // namespace kynema::solver
Definition calculate_error_sum_squares.hpp:5
A Kernel which sums the system matrix contributions computed at a mass element's node into the correc...
Definition contribute_masses_to_sparse_matrix.hpp:15
CrsMatrixType sparse
Definition contribute_masses_to_sparse_matrix.hpp:30
typename CrsMatrixType::device_type DeviceType
Definition contribute_masses_to_sparse_matrix.hpp:16
typename TeamPolicy::member_type member_type
Definition contribute_masses_to_sparse_matrix.hpp:20
ConstView< size_t *[6]> element_freedom_table
Definition contribute_masses_to_sparse_matrix.hpp:28
double conditioner
Definition contribute_masses_to_sparse_matrix.hpp:26
ConstView< dof::FreedomSignature * > element_freedom_signature
Definition contribute_masses_to_sparse_matrix.hpp:27
typename View< ValueType >::const_type ConstView
Definition contribute_masses_to_sparse_matrix.hpp:24
typename CrsMatrixType::staticcrsgraph_type::entries_type::non_const_type ColIdxType
Definition contribute_masses_to_sparse_matrix.hpp:18
KOKKOS_FUNCTION void operator()(member_type member) const
Definition contribute_masses_to_sparse_matrix.hpp:33
typename CrsMatrixType::values_type::non_const_type RowDataType
Definition contribute_masses_to_sparse_matrix.hpp:17
Kokkos::TeamPolicy< typename DeviceType::execution_space > TeamPolicy
Definition contribute_masses_to_sparse_matrix.hpp:19
ConstView< double *[6][6]> dense
Definition contribute_masses_to_sparse_matrix.hpp:29
Kokkos::View< ValueType, DeviceType > View
Definition contribute_masses_to_sparse_matrix.hpp:22