Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_getGraphDiagOffsets_def.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_DEF_HPP
45 #define TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DEF_HPP
46 
51 
53 #include "Tpetra_Map.hpp"
54 #include "KokkosSparse_findRelOffset.hpp"
55 
56 namespace Tpetra {
57 namespace Details {
58 namespace Impl {
59 
73 template<class LO,
74  class GO,
75  class Node,
76  class DiagOffsetType>
78 GetGraphDiagOffsets (const diag_offsets_type& diagOffsets,
79  const local_map_type& lclRowMap,
80  const local_map_type& lclColMap,
81  const row_offsets_type& ptr,
82  const lcl_col_inds_type& ind,
83  const bool isSorted) :
84  diagOffsets_ (diagOffsets),
85  lclRowMap_ (lclRowMap),
86  lclColMap_ (lclColMap),
87  ptr_ (ptr),
88  ind_ (ind),
89  isSorted_ (isSorted)
90 {
91  typedef typename device_type::execution_space execution_space;
92  typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
93 
94  const LO lclNumRows = lclRowMap.getNodeNumElements ();
95  policy_type range (0, lclNumRows);
96  Kokkos::parallel_for (range, *this);
97 }
98 
99 template<class LO,
100  class GO,
101  class Node,
102  class DiagOffsetType>
103 KOKKOS_FUNCTION void
105 operator() (const LO& lclRowInd) const
106 {
107  const size_t STINV =
108  Tpetra::Details::OrdinalTraits<diag_offset_type>::invalid ();
109  const GO gblRowInd = lclRowMap_.getGlobalElement (lclRowInd);
110  const GO gblColInd = gblRowInd;
111  const LO lclColInd = lclColMap_.getLocalElement (gblColInd);
112 
113  if (lclColInd == Tpetra::Details::OrdinalTraits<LO>::invalid ()) {
114  diagOffsets_[lclRowInd] = STINV;
115  }
116  else {
117  // Could be empty, but that's OK.
118  const LO numEnt = ptr_[lclRowInd+1] - ptr_[lclRowInd];
119  // std::pair doesn't have its methods marked as device
120  // functions, so we have to use Kokkos::pair.
121  auto lclColInds =
122  Kokkos::subview (ind_, Kokkos::make_pair (ptr_[lclRowInd],
123  ptr_[lclRowInd+1]));
124  using ::KokkosSparse::findRelOffset;
125  const LO diagOffset =
126  findRelOffset<LO, lcl_col_inds_type> (lclColInds, numEnt,
127  lclColInd, 0, isSorted_);
128  diagOffsets_[lclRowInd] = (diagOffset == numEnt) ? STINV :
129  static_cast<diag_offset_type> (diagOffset);
130  }
131 }
132 
133 } // namespace Impl
134 } // namespace Details
135 } // namespace Tpetra
136 
137 // Explicit template instantiation macro for
138 // Tpetra::Details::Impl::GetGraphDiagOffsets. NOT FOR USERS!!! Must
139 // be used inside the Tpetra namespace.
140 #define TPETRA_DETAILS_IMPL_GETGRAPHDIAGOFFSETS_INSTANT( LO, GO, NODE ) \
141  template class Details::Impl::GetGraphDiagOffsets< LO, GO, NODE::device_type >;
142 
143 #endif // TPETRA_DETAILS_GETGRAPHDIAGOFFSETS_DEF_HPP
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.
Import KokkosSparse::OrdinalTraits, a traits class for &quot;invalid&quot; (flag) values of integer types...
KOKKOS_INLINE_FUNCTION LocalOrdinal getNodeNumElements() const
The number of indices that live on the calling process.
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.