Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_getEntryOnHost.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_GETENTRYONHOST_HPP
11#define TPETRA_DETAILS_GETENTRYONHOST_HPP
12
18
19#include "TpetraCore_config.h"
20#include "Kokkos_Core.hpp"
21
22namespace Tpetra {
23namespace Details {
24
25template<class ViewType,
26 class IndexType>
27typename ViewType::non_const_value_type
28getEntryOnHost (const ViewType& x,
29 const IndexType ind)
30{
31 using execution_space = typename ViewType::execution_space;
32 static_assert (ViewType::rank == 1, "x must be a rank-1 Kokkos::View.");
33 // Get a 0-D subview of the entry of the array, and copy to host scalar.
34 typename ViewType::non_const_value_type val;
35 // DEEP_COPY REVIEW - DEVICE-TO-VALUE
36 Kokkos::deep_copy(execution_space(), val, Kokkos::subview(x, ind));
37 return val;
38}
39
40template<class ViewType,
41 class IndexType>
42typename ViewType::HostMirror::const_type
43getEntriesOnHost (const ViewType& x,
44 const IndexType ind,
45 const int count)
46{
47 static_assert (ViewType::rank == 1, "x must be a rank-1 Kokkos::View.");
48 // Get a 1-D subview of the entry of the array, and copy to host.
49 auto subview = Kokkos::subview(x, Kokkos::make_pair(ind, ind + count));
50 return Kokkos::create_mirror_view_and_copy(typename ViewType::HostMirror::memory_space(), subview);
51}
52
53} // namespace Details
54} // namespace Tpetra
55
56#endif // TPETRA_DETAILS_GETENTRYONHOST_HPP
Nonmember function that computes a residual Computes R = B - A * X.
Namespace Tpetra contains the class and methods constituting the Tpetra library.