Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_getGraphDiagOffsets_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DECL_HPP
11#define TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DECL_HPP
12
17
18#include "TpetraCore_config.h"
19#include "Kokkos_Core.hpp"
20#include "KokkosSparse_StaticCrsGraph.hpp"
22#include <type_traits>
23
24namespace Tpetra {
25namespace Details {
26namespace Impl {
27
41template <class LO,
42 class GO,
43 class DeviceType,
44 class DiagOffsetType = size_t>
46 public:
47 typedef typename DeviceType::device_type device_type;
49 typedef ::Kokkos::View<diag_offset_type*,
50 device_type,
51 ::Kokkos::MemoryUnmanaged>
52 diag_offsets_type;
53 using local_graph_device_type = ::KokkosSparse::StaticCrsGraph<LO,
54 ::Kokkos::LayoutLeft,
55 device_type, void, size_t>;
56 typedef ::Tpetra::Details::LocalMap<LO, GO, device_type> local_map_type;
57 typedef ::Kokkos::View<const typename local_graph_device_type::size_type*,
58 ::Kokkos::LayoutLeft,
59 device_type,
60 ::Kokkos::MemoryUnmanaged>
61 row_offsets_type;
62 // This is unmanaged for performance, because we need to take
63 // subviews inside the functor.
64 typedef ::Kokkos::View<const LO*,
65 ::Kokkos::LayoutLeft,
66 device_type,
67 ::Kokkos::MemoryUnmanaged>
68 lcl_col_inds_type;
69
71 GetGraphDiagOffsets(const diag_offsets_type& diagOffsets,
74 const row_offsets_type& ptr,
75 const lcl_col_inds_type& ind,
76 const bool isSorted);
77
79 KOKKOS_FUNCTION void operator()(const LO& lclRowInd) const;
80
81 private:
82 diag_offsets_type diagOffsets_;
83 local_map_type lclRowMap_;
84 local_map_type lclColMap_;
85 row_offsets_type ptr_;
86 lcl_col_inds_type ind_;
87 bool isSorted_;
88};
89
90} // namespace Impl
91
92template <class DiagOffsetsType,
93 class LclMapType,
94 class RowOffsetsType,
95 class LclColIndsType>
96void getGraphDiagOffsets(const DiagOffsetsType& diagOffsets,
97 const LclMapType& lclRowMap,
98 const LclMapType& lclColMap,
99 const RowOffsetsType& ptr,
100 const LclColIndsType& ind,
101 const bool isSorted) {
102 static_assert(Kokkos::is_view<DiagOffsetsType>::value,
103 "DiagOffsetsType (the type of diagOffsets) must be a Kokkos::View.");
104 static_assert(Kokkos::is_view<RowOffsetsType>::value,
105 "RowOffsetsType (the type of ptr) must be a Kokkos::View.");
106 static_assert(Kokkos::is_view<LclColIndsType>::value,
107 "LclColIndsType (the type of ind) must be a Kokkos::View.");
108 static_assert(static_cast<int>(DiagOffsetsType::rank) == 1,
109 "DiagOffsetsType (the type of diagOffsets) must be a rank-1 Kokkos::View.");
110 static_assert(static_cast<int>(RowOffsetsType::rank) == 1,
111 "RowOffsetsType (the type of ptr) must be a rank-1 Kokkos::View.");
112 static_assert(static_cast<int>(LclColIndsType::rank) == 1,
113 "LclColIndsType (the type of ind) must be a rank-1 Kokkos::View.");
114 typedef typename DiagOffsetsType::non_const_value_type diag_offset_type;
115 static_assert(std::is_same<typename DiagOffsetsType::value_type, diag_offset_type>::value,
116 "DiagOffsetsType (the type of diagOffsets) must be a nonconst "
117 "Kokkos::View, since it is the output argument of this function.");
118 static_assert(std::is_integral<diag_offset_type>::value,
119 "The type of each entry of diagOffsets must be an integer.");
120 typedef typename LclColIndsType::non_const_value_type local_ordinal_type;
121 static_assert(std::is_integral<local_ordinal_type>::value,
122 "The type of each entry of ind (the array of column indices) "
123 "must be an integer.");
124 static_assert(sizeof(diag_offset_type) >= sizeof(local_ordinal_type),
125 "Diagonal offset type must be big enough to count the number "
126 "of column indices, since the diagonal entry (if it exists) "
127 "may be anywhere in a row.");
128 typedef typename RowOffsetsType::non_const_value_type row_offset_type;
129 static_assert(std::is_integral<row_offset_type>::value,
130 "The type of each entry of ptr must be an integer.");
131
132 typedef typename LclMapType::local_ordinal_type LO;
133 typedef typename LclMapType::global_ordinal_type GO;
134 typedef typename LclMapType::device_type DT;
135
137 // The functor's constructor runs the functor.
138 impl_type impl(diagOffsets, lclRowMap, lclColMap, ptr, ind, isSorted);
139}
140
141} // namespace Details
142} // namespace Tpetra
143
144#endif // TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DECL_HPP
Declaration and definition of the Tpetra::Map class, an implementation detail of Tpetra::Map.
Struct that holds views of the contents of a CrsMatrix.
Implementation detail of Tpetra::Details::getGraphDiagOffsets, which in turn is an implementation det...
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.
Implementation details of Tpetra.
Namespace Tpetra contains the class and methods constituting the Tpetra library.