42#ifndef TPETRA_DETAILS_FILL_HPP
43#define TPETRA_DETAILS_FILL_HPP
67memsetWrapper (
void* dest,
int ch, std::size_t count);
70template<
class ViewType,
74 const int rank = ViewType::Rank>
84 static_assert (std::is_integral<IndexType>::value,
85 "IndexType must be a built-in integer type.");
100 static_assert (ViewType::Rank == 1,
101 "ViewType must be a rank-1 Kokkos::View.");
102 static_assert (std::is_integral<IndexType>::value,
103 "IndexType must be a built-in integer type.");
119 typedef Kokkos::RangePolicy<ExecutionSpace, IndexType> range_type;
120 Kokkos::parallel_for (
"fill", range_type (0,
numRows),
139 X_ (
X), alpha_ (
alpha), numCols_ (numCols)
141 static_assert (ViewType::Rank == 2,
142 "ViewType must be a rank-2 Kokkos::View.");
143 static_assert (std::is_integral<IndexType>::value,
144 "IndexType must be a built-in integer type.");
162 typedef Kokkos::RangePolicy<ExecutionSpace, IndexType> range_type;
163 Kokkos::parallel_for (
"fill", range_type (0,
numRows),
173#if defined(KOKKOS_ENABLE_SERIAL)
186 fill (
const Kokkos::Serial& ,
192 static_assert (ViewType::Rank == 1,
193 "ViewType must be a rank-1 Kokkos::View.");
194 static_assert (std::is_integral<IndexType>::value,
195 "IndexType must be a built-in integer type.");
196 using ::Tpetra::Details::Blas::BlasSupportsScalar;
215#if defined(KOKKOS_ENABLE_SERIAL)
218template<
class ViewType,
228 fill (
const Kokkos::Serial& ,
230 const ValueType& alpha,
232 const IndexType numCols)
234 static_assert (ViewType::Rank == 2,
235 "ViewType must be a rank-2 Kokkos::View.");
236 static_assert (std::is_integral<IndexType>::value,
237 "IndexType must be a built-in integer type.");
238 using ::Tpetra::Details::Blas::BlasSupportsScalar;
239 typedef typename ViewType::non_const_value_type view_value_type;
240 typedef typename ViewType::array_layout array_layout;
244 constexpr bool podType = BlasSupportsScalar<view_value_type>::value;
246 if (podType && alpha == ValueType (0.0)) {
247 if (X.span_is_contiguous ()) {
248 memsetWrapper (X.data (), 0, X.span () * sizeof (view_value_type));
250 else if (std::is_same<array_layout, Kokkos::LayoutLeft>::value) {
252 for (IndexType j = 0; j < numCols; ++j) {
253 auto X_j = Kokkos::subview (X, Kokkos::ALL (), j);
254 memsetWrapper (X_j.data (), 0,
255 X_j.extent (0) * sizeof (view_value_type));
259 Kokkos::deep_copy (X, view_value_type (0.0));
263 Kokkos::deep_copy (X, alpha);
282template<
class ViewType,
285 class ExecutionSpace>
293 static_assert (std::is_integral<IndexType>::value,
294 "IndexType must be a built-in integer type.");
312 static_assert (ViewType::Rank == 2,
"ViewType must be a rank-2 "
313 "Kokkos::View in order to call the \"whichVectors\" "
314 "specialization of fill.");
315 static_assert (std::is_integral<IndexType>::value,
316 "IndexType must be a built-in integer type.");
319 auto X_j = Kokkos::subview (
X, Kokkos::ALL (),
j);
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 of Tpetra::Details::Blas::fill.
Implementation details of Tpetra.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Do BLAS libraries (all that are compliant with the BLAS Standard) support the given "scalar" (matrix ...