Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
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#if KOKKOS_VERSION >= 40799
25#include "KokkosKernels_ArithTraits.hpp"
26#else
27#include "Kokkos_ArithTraits.hpp"
28#endif
30#include "Tpetra_RowMatrix_decl.hpp"
32#include <type_traits>
33
34namespace Tpetra {
35namespace Details {
36
49template <class DiagType,
50 class LocalMapType,
51 class CrsMatrixType>
53 typedef typename LocalMapType::local_ordinal_type LO; // local ordinal type
54 typedef typename LocalMapType::global_ordinal_type GO; // global ordinal type
55 typedef typename CrsMatrixType::device_type device_type;
56 typedef typename CrsMatrixType::value_type scalar_type;
57 typedef typename CrsMatrixType::size_type offset_type;
58
60 typedef LO value_type;
61
70 const LocalMapType& rowMap,
71 const LocalMapType& colMap,
72 const CrsMatrixType& A)
73 : D_(D)
74 , rowMap_(rowMap)
75 , colMap_(colMap)
76 , A_(A) {}
77
83 const LO INV = Tpetra::Details::OrdinalTraits<LO>::invalid();
84 const scalar_type ZERO =
85#if KOKKOS_VERSION >= 40799
86 KokkosKernels::ArithTraits<scalar_type>::zero();
87#else
88 Kokkos::ArithTraits<scalar_type>::zero();
89#endif
90
91 // If the row lacks a stored diagonal entry, then its value is zero.
92 D_(lclRowInd) = ZERO;
93 const GO gblInd = rowMap_.getGlobalElement(lclRowInd);
94 const LO lclColInd = colMap_.getLocalElement(gblInd);
95
96 if (lclColInd != INV) {
97 auto curRow = A_.rowConst(lclRowInd);
98
99 // FIXME (mfh 12 May 2016) Use binary search when the row is
100 // long enough. findRelOffset currently lives in KokkosKernels
101 // (in tpetra/kernels/src/Kokkos_Sparse_findRelOffset.hpp).
102 LO offset = 0;
103 const LO numEnt = curRow.length;
104 for (; offset < numEnt; ++offset) {
105 if (curRow.colidx(offset) == lclColInd) {
106 break;
107 }
108 }
109
110 if (offset == numEnt) {
111 ++errCount;
112 } else {
113 D_(lclRowInd) = curRow.value(offset);
114 }
115 }
116 }
117
118 private:
120 DiagType D_;
122 LocalMapType rowMap_;
124 LocalMapType colMap_;
126 CrsMatrixType A_;
127};
128
151template <class DiagType,
152 class LocalMapType,
153 class CrsMatrixType>
154static typename LocalMapType::local_ordinal_type
156 const LocalMapType& rowMap,
157 const LocalMapType& colMap,
158 const CrsMatrixType& A) {
159 static_assert(Kokkos::is_view<DiagType>::value,
160 "DiagType must be a Kokkos::View.");
161 static_assert(static_cast<int>(DiagType::rank) == 1,
162 "DiagType must be a 1-D Kokkos::View.");
163 static_assert(std::is_same<typename DiagType::value_type, typename DiagType::non_const_value_type>::value,
164 "DiagType must be a nonconst Kokkos::View.");
165 typedef typename LocalMapType::local_ordinal_type LO;
166 typedef typename CrsMatrixType::device_type device_type;
167 typedef typename device_type::execution_space execution_space;
168 typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
169
170 typedef Kokkos::View<typename DiagType::non_const_value_type*,
171 typename DiagType::array_layout,
172 typename DiagType::device_type,
173 Kokkos::MemoryUnmanaged>
174 diag_type;
175 diag_type D_out = D;
177 functor(D_out, rowMap, colMap, A);
178 const LO numRows = static_cast<LO>(D.extent(0));
179 LO errCount = 0;
180 Kokkos::parallel_reduce(policy_type(0, numRows), functor, errCount);
181 return errCount;
182}
183
220template <class SC, class LO, class GO, class NT>
222 const ::Tpetra::RowMatrix<SC, LO, GO, NT>& A,
223 const bool debug =
225 true);
226#else // ! HAVE_TPETRA_DEBUG
227 false);
228#endif // HAVE_TPETRA_DEBUG
229
230} // namespace Details
231} // namespace Tpetra
232
233#endif // TPETRA_DETAILS_GETDIAGCOPYWITHOUTOFFSETS_DECL_HPP
Import KokkosSparse::OrdinalTraits, a traits class for "invalid" (flag) values of integer types,...
Declaration of the Tpetra::Vector class.
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
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's diagonal entries into a Tpetra::V...
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,...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
@ ZERO
Replace old values with zero.
Functor that implements much of the one-argument overload of Tpetra::CrsMatrix::getLocalDiagCopy,...
CrsMatrixGetDiagCopyFunctor(const DiagType &D, const LocalMapType &rowMap, const LocalMapType &colMap, const CrsMatrixType &A)
Constructor.
KOKKOS_FUNCTION void operator()(const LO &lclRowInd, value_type &errCount) const
Operator for Kokkos::parallel_for.