Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_fill.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_FILL_HPP
11#define TPETRA_DETAILS_FILL_HPP
12
21
23#include <type_traits>
24
25namespace Tpetra {
26namespace Details {
27namespace Blas {
28
36template <class ViewType,
37 class ValueType,
38 class IndexType,
39 class ExecutionSpace>
41 const ViewType& X,
42 const ValueType& alpha,
43 const IndexType numRows,
44 const IndexType numCols) {
45 static_assert(std::is_integral<IndexType>::value,
46 "IndexType must be a built-in integer type.");
47 auto X_j = Kokkos::subview(X, Kokkos::make_pair(IndexType(0), numRows), Kokkos::make_pair(IndexType(0), numCols));
48 Kokkos::deep_copy(execSpace, X_j, alpha);
49}
50
51template <class ViewType,
52 class ValueType,
53 class IndexType,
54 class ExecutionSpace>
55void fill(const ExecutionSpace& execSpace,
56 const ViewType& X,
57 const ValueType& alpha,
58 const IndexType numRows,
59 const IndexType numCols,
60 const size_t whichVectors[]) {
61 static_assert(ViewType::rank == 2,
62 "ViewType must be a rank-2 "
63 "Kokkos::View in order to call the \"whichVectors\" "
64 "specialization of fill.");
65 static_assert(std::is_integral<IndexType>::value,
66 "IndexType must be a built-in integer type.");
67 for (IndexType k = 0; k < numCols; ++k) {
68 const IndexType j = whichVectors[k];
69 auto X_j = Kokkos::subview(X, Kokkos::make_pair(IndexType(0), numRows), j);
70 Kokkos::deep_copy(execSpace, X_j, alpha);
71 }
72}
73
74} // namespace Blas
75} // namespace Details
76} // namespace Tpetra
77
78#endif // TPETRA_DETAILS_FILL_HPP
Type traits for Tpetra's BLAS wrappers; an implementation detail of Tpetra::MultiVector.
void fill(const ExecutionSpace &execSpace, const ViewType &X, const ValueType &alpha, const IndexType numRows, const IndexType numCols)
Fill the entries of the given 1-D or 2-D Kokkos::View with the given scalar value alpha.
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
Namespace Tpetra contains the class and methods constituting the Tpetra library.