Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_getGraphDiagOffsets_decl.hpp
Go to the documentation of this file.
1 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // Tpetra: Templated Linear Algebra Services Package
6 // Copyright (2008) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 */
43 
44 #ifndef TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DECL_HPP
45 #define TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DECL_HPP
46 
51 
52 #include "TpetraCore_config.h"
53 #include "Kokkos_Core.hpp"
54 #include "Kokkos_StaticCrsGraph.hpp"
56 #include <type_traits>
57 
58 namespace Tpetra {
59 namespace Details {
60 namespace Impl {
61 
75 template<class LO,
76  class GO,
77  class DeviceType,
78  class DiagOffsetType = size_t>
80 public:
81  typedef typename DeviceType::device_type device_type;
82  typedef DiagOffsetType diag_offset_type;
83  typedef ::Kokkos::View<diag_offset_type*,
84  device_type,
85  ::Kokkos::MemoryUnmanaged> diag_offsets_type;
86  typedef ::Kokkos::StaticCrsGraph<LO,
87  ::Kokkos::LayoutLeft,
88  device_type> local_graph_type;
89  typedef ::Tpetra::Details::LocalMap<LO, GO, device_type> local_map_type;
90  typedef ::Kokkos::View<const typename local_graph_type::size_type*,
91  ::Kokkos::LayoutLeft,
92  device_type,
93  ::Kokkos::MemoryUnmanaged> row_offsets_type;
94  // This is unmanaged for performance, because we need to take
95  // subviews inside the functor.
96  typedef ::Kokkos::View<const LO*,
97  ::Kokkos::LayoutLeft,
98  device_type,
99  ::Kokkos::MemoryUnmanaged> lcl_col_inds_type;
100 
102  GetGraphDiagOffsets (const diag_offsets_type& diagOffsets,
103  const local_map_type& lclRowMap,
104  const local_map_type& lclColMap,
105  const row_offsets_type& ptr,
106  const lcl_col_inds_type& ind,
107  const bool isSorted);
108 
110  KOKKOS_FUNCTION void operator() (const LO& lclRowInd) const;
111 
112 private:
113  diag_offsets_type diagOffsets_;
114  local_map_type lclRowMap_;
115  local_map_type lclColMap_;
116  row_offsets_type ptr_;
117  lcl_col_inds_type ind_;
118  bool isSorted_;
119 };
120 
121 } // namespace Impl
122 
123 template<class DiagOffsetsType,
124  class LclMapType,
125  class RowOffsetsType,
126  class LclColIndsType>
127 void
128 getGraphDiagOffsets (const DiagOffsetsType& diagOffsets,
129  const LclMapType& lclRowMap,
130  const LclMapType& lclColMap,
131  const RowOffsetsType& ptr,
132  const LclColIndsType& ind,
133  const bool isSorted)
134 {
135  static_assert (Kokkos::Impl::is_view<DiagOffsetsType>::value,
136  "DiagOffsetsType (the type of diagOffsets) must be a Kokkos::View.");
137  static_assert (Kokkos::Impl::is_view<RowOffsetsType>::value,
138  "RowOffsetsType (the type of ptr) must be a Kokkos::View.");
139  static_assert (Kokkos::Impl::is_view<LclColIndsType>::value,
140  "LclColIndsType (the type of ind) must be a Kokkos::View.");
141  static_assert (static_cast<int> (DiagOffsetsType::rank) == 1,
142  "DiagOffsetsType (the type of diagOffsets) must be a rank-1 Kokkos::View.");
143  static_assert (static_cast<int> (RowOffsetsType::rank) == 1,
144  "RowOffsetsType (the type of ptr) must be a rank-1 Kokkos::View.");
145  static_assert (static_cast<int> (LclColIndsType::rank) == 1,
146  "LclColIndsType (the type of ind) must be a rank-1 Kokkos::View.");
147  typedef typename DiagOffsetsType::non_const_value_type diag_offset_type;
148  static_assert (std::is_same<typename DiagOffsetsType::value_type, diag_offset_type>::value,
149  "DiagOffsetsType (the type of diagOffsets) must be a nonconst "
150  "Kokkos::View, since it is the output argument of this function.");
151  static_assert (std::is_integral<diag_offset_type>::value,
152  "The type of each entry of diagOffsets must be an integer.");
153  typedef typename LclColIndsType::non_const_value_type local_ordinal_type;
154  static_assert (std::is_integral<local_ordinal_type>::value,
155  "The type of each entry of ind (the array of column indices) "
156  "must be an integer.");
157  static_assert (sizeof (diag_offset_type) >= sizeof (local_ordinal_type),
158  "Diagonal offset type must be big enough to count the number "
159  "of column indices, since the diagonal entry (if it exists) "
160  "may be anywhere in a row.");
161  typedef typename RowOffsetsType::non_const_value_type row_offset_type;
162  static_assert (std::is_integral<row_offset_type>::value,
163  "The type of each entry of ptr must be an integer.");
164 
165  typedef typename LclMapType::local_ordinal_type LO;
166  typedef typename LclMapType::global_ordinal_type GO;
167  typedef typename LclMapType::device_type DT;
168 
170  // The functor's constructor runs the functor.
171  impl_type impl (diagOffsets, lclRowMap, lclColMap, ptr, ind, isSorted);
172 }
173 
174 } // namespace Details
175 } // namespace Tpetra
176 
177 #endif // TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DECL_HPP
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.
Declaration and definition of the Tpetra::Map class, an implementation detail of Tpetra::Map.
Implementation detail of Tpetra::Details::getGraphDiagOffsets, which in turn is an implementation det...
GetGraphDiagOffsets(const diag_offsets_type &diagOffsets, const local_map_type &lclRowMap, const local_map_type &lclColMap, const row_offsets_type &ptr, const lcl_col_inds_type &ind, const bool isSorted)
Constructor; also runs the functor.