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#include "KokkosKernels_ArithTraits.hpp"
26#include "Tpetra_RowMatrix_decl.hpp"
28#include <type_traits>
29
30namespace Tpetra {
31namespace Details {
32
45template <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
66 const LocalMapType& rowMap,
67 const LocalMapType& colMap,
68 const CrsMatrixType& A)
69 : D_(D)
70 , rowMap_(rowMap)
71 , colMap_(colMap)
72 , A_(A) {}
73
80 const LO INV = Tpetra::Details::OrdinalTraits<LO>::invalid();
81 const scalar_type ZERO =
82 KokkosKernels::ArithTraits<scalar_type>::zero();
83
84 // If the row lacks a stored diagonal entry, then its value is zero.
85 D_(lclRowInd) = ZERO;
86 const GO gblInd = rowMap_.getGlobalElement(lclRowInd);
87 const LO lclColInd = colMap_.getLocalElement(gblInd);
88
89 if (lclColInd != INV) {
90 auto curRow = A_.rowConst(lclRowInd);
91
92 // FIXME (mfh 12 May 2016) Use binary search when the row is
93 // long enough. findRelOffset currently lives in KokkosKernels
94 // (in tpetra/kernels/src/Kokkos_Sparse_findRelOffset.hpp).
95 LO offset = 0;
96 const LO numEnt = curRow.length;
97 for (; offset < numEnt; ++offset) {
98 if (curRow.colidx(offset) == lclColInd) {
99 break;
100 }
101 }
102
103 if (offset == numEnt) {
104 ++errCount;
105 } else {
106 D_(lclRowInd) = curRow.value(offset);
107 }
108 }
109 }
110
111 private:
113 DiagType D_;
115 LocalMapType rowMap_;
117 LocalMapType colMap_;
119 CrsMatrixType A_;
120};
121
144template <class DiagType,
145 class LocalMapType,
146 class CrsMatrixType>
147static typename LocalMapType::local_ordinal_type
149 const LocalMapType& rowMap,
150 const LocalMapType& colMap,
151 const CrsMatrixType& A) {
152 static_assert(Kokkos::is_view<DiagType>::value,
153 "DiagType must be a Kokkos::View.");
154 static_assert(static_cast<int>(DiagType::rank) == 1,
155 "DiagType must be a 1-D Kokkos::View.");
156 static_assert(std::is_same<typename DiagType::value_type, typename DiagType::non_const_value_type>::value,
157 "DiagType must be a nonconst Kokkos::View.");
158 typedef typename LocalMapType::local_ordinal_type LO;
159 typedef typename CrsMatrixType::device_type device_type;
160 typedef typename device_type::execution_space execution_space;
161 typedef Kokkos::RangePolicy<execution_space, LO> policy_type;
162
163 typedef Kokkos::View<typename DiagType::non_const_value_type*,
164 typename DiagType::array_layout,
165 typename DiagType::device_type,
166 Kokkos::MemoryUnmanaged>
167 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
213template <class SC, class LO, class GO, class NT>
215 const ::Tpetra::RowMatrix<SC, LO, GO, NT>& A,
216 const bool debug =
218 true);
219#else // ! HAVE_TPETRA_DEBUG
220 false);
221#endif // HAVE_TPETRA_DEBUG
222
223} // namespace Details
224} // namespace Tpetra
225
226#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.