Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_getDiagCopyWithoutOffsets_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_GETDIAGCOPYWITHOUTOFFSETS_DECL_HPP
11 #define TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DECL_HPP
12 
21 
22 #include "TpetraCore_config.h"
23 #include "Kokkos_Core.hpp"
24 #include "Kokkos_ArithTraits.hpp"
26 #include "Tpetra_RowMatrix_decl.hpp"
27 #include "Tpetra_Vector_decl.hpp"
28 #include <type_traits>
29 
30 namespace Tpetra {
31 namespace Details {
32 
45 template<class DiagType,
46  class LocalMapType,
47  class CrsMatrixType>
49  typedef typename LocalMapType::local_ordinal_type LO; // local ordinal type
50  typedef typename LocalMapType::global_ordinal_type GO; // global ordinal type
51  typedef typename CrsMatrixType::device_type device_type;
52  typedef typename CrsMatrixType::value_type scalar_type;
53  typedef typename CrsMatrixType::size_type offset_type;
54 
56  typedef LO value_type;
57 
65  CrsMatrixGetDiagCopyFunctor (const DiagType& D,
66  const LocalMapType& rowMap,
67  const LocalMapType& colMap,
68  const CrsMatrixType& A) :
69  D_ (D), rowMap_ (rowMap), colMap_ (colMap), A_ (A)
70  {}
71 
75  KOKKOS_FUNCTION void
76  operator () (const LO& lclRowInd, value_type& errCount) const
77  {
78  const LO INV = Tpetra::Details::OrdinalTraits<LO>::invalid ();
79  const scalar_type ZERO =
80  Kokkos::ArithTraits<scalar_type>::zero ();
81 
82  // If the row lacks a stored diagonal entry, then its value is zero.
83  D_(lclRowInd) = ZERO;
84  const GO gblInd = rowMap_.getGlobalElement (lclRowInd);
85  const LO lclColInd = colMap_.getLocalElement (gblInd);
86 
87  if (lclColInd != INV) {
88  auto curRow = A_.rowConst (lclRowInd);
89 
90  // FIXME (mfh 12 May 2016) Use binary search when the row is
91  // long enough. findRelOffset currently lives in KokkosKernels
92  // (in tpetra/kernels/src/Kokkos_Sparse_findRelOffset.hpp).
93  LO offset = 0;
94  const LO numEnt = curRow.length;
95  for ( ; offset < numEnt; ++offset) {
96  if (curRow.colidx(offset) == lclColInd) {
97  break;
98  }
99  }
100 
101  if (offset == numEnt) {
102  ++errCount;
103  }
104  else {
105  D_(lclRowInd) = curRow.value(offset);
106  }
107  }
108  }
109 
110 private:
112  DiagType D_;
114  LocalMapType rowMap_;
116  LocalMapType colMap_;
118  CrsMatrixType A_;
119 };
120 
121 
144 template<class DiagType,
145  class LocalMapType,
146  class CrsMatrixType>
147 static typename LocalMapType::local_ordinal_type
148 getDiagCopyWithoutOffsets (const DiagType& D,
149  const LocalMapType& rowMap,
150  const LocalMapType& colMap,
151  const CrsMatrixType& A)
152 {
153  static_assert (Kokkos::is_view<DiagType>::value,
154  "DiagType must be a Kokkos::View.");
155  static_assert (static_cast<int> (DiagType::rank) == 1,
156  "DiagType must be a 1-D Kokkos::View.");
157  static_assert (std::is_same<typename DiagType::value_type, typename DiagType::non_const_value_type>::value,
158  "DiagType must be a nonconst Kokkos::View.");
159  typedef typename LocalMapType::local_ordinal_type LO;
160  typedef typename CrsMatrixType::device_type device_type;
161  typedef typename device_type::execution_space execution_space;
162  typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
163 
164  typedef Kokkos::View<typename DiagType::non_const_value_type*,
165  typename DiagType::array_layout,
166  typename DiagType::device_type,
167  Kokkos::MemoryUnmanaged> diag_type;
168  diag_type D_out = D;
170  functor (D_out, rowMap, colMap, A);
171  const LO numRows = static_cast<LO> (D.extent (0));
172  LO errCount = 0;
173  Kokkos::parallel_reduce (policy_type (0, numRows), functor, errCount);
174  return errCount;
175 }
176 
213 template<class SC, class LO, class GO, class NT>
214 LO
216  const ::Tpetra::RowMatrix<SC, LO, GO, NT>& A,
217  const bool debug =
218 #ifdef HAVE_TPETRA_DEBUG
219  true);
220 #else // ! HAVE_TPETRA_DEBUG
221  false);
222 #endif // HAVE_TPETRA_DEBUG
223 
224 } // namespace Details
225 } // namespace Tpetra
226 
227 #endif // TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DECL_HPP
Import KokkosSparse::OrdinalTraits, a traits class for &quot;invalid&quot; (flag) values of integer types...
Functor that implements much of the one-argument overload of Tpetra::CrsMatrix::getLocalDiagCopy, for the case where the matrix is fill complete.
LO value_type
The result of the reduction; number of errors.
Declaration of the Tpetra::Vector class.
KOKKOS_FUNCTION void operator()(const LO &lclRowInd, value_type &errCount) const
Operator for Kokkos::parallel_for.
static LocalMapType::local_ordinal_type getDiagCopyWithoutOffsets(const DiagType &D, const LocalMapType &rowMap, const LocalMapType &colMap, const CrsMatrixType &A)
Given a locally indexed, local sparse matrix, and corresponding local row and column Maps...
LO getLocalDiagCopyWithoutOffsetsNotFillComplete(::Tpetra::Vector< SC, LO, GO, NT > &diag, const ::Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool debug=false)
Given a locally indexed, global sparse matrix, extract the matrix&#39;s diagonal entries into a Tpetra::V...
Replace old values with zero.
A distributed dense vector.
CrsMatrixGetDiagCopyFunctor(const DiagType &D, const LocalMapType &rowMap, const LocalMapType &colMap, const CrsMatrixType &A)
Constructor.