42 #ifndef TPETRA_DETAILS_FILL_HPP 
   43 #define TPETRA_DETAILS_FILL_HPP 
   58 #include <type_traits> 
   67 memsetWrapper (
void* dest, 
int ch, std::size_t count);
 
   70 template<
class ViewType,
 
   74          const int rank = ViewType::Rank>
 
   78   fill  (
const ExecutionSpace& execSpace,
 
   80          const ValueType& alpha,
 
   84     static_assert (std::is_integral<IndexType>::value,
 
   85                    "IndexType must be a built-in integer type.");
 
   91 template<
class ViewType,
 
   95 class Fill<ViewType, ValueType, ExecutionSpace, IndexType, 1> {
 
   97   Fill (
const ViewType& X, 
const ValueType& alpha) :
 
   98     X_ (X), alpha_ (alpha)
 
  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.");
 
  106   KOKKOS_INLINE_FUNCTION 
void 
  107   operator () (
const IndexType& i)
 const 
  113   fill (
const ExecutionSpace& execSpace,
 
  115         const ValueType& alpha,
 
  116         const IndexType numRows,
 
  119     typedef Kokkos::RangePolicy<ExecutionSpace, IndexType> range_type;
 
  120     Kokkos::parallel_for (
"fill", range_type (0, numRows),
 
  130 template<
class ViewType,
 
  132          class ExecutionSpace,
 
  134 class Fill<ViewType, ValueType, ExecutionSpace, IndexType, 2> {
 
  136   Fill (
const ViewType& X,
 
  137         const ValueType& alpha,
 
  138         const IndexType numCols) :
 
  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.");
 
  147   KOKKOS_INLINE_FUNCTION 
void 
  148   operator () (
const IndexType& i)
 const 
  150     for (IndexType j = 0; j < numCols_; ++j) {
 
  156   fill (
const ExecutionSpace& execSpace,
 
  158         const ValueType& alpha,
 
  159         const IndexType numRows,
 
  160         const IndexType numCols)
 
  162     typedef Kokkos::RangePolicy<ExecutionSpace, IndexType> range_type;
 
  163     Kokkos::parallel_for (
"fill", range_type (0, numRows),
 
  164                           Fill (X, alpha, numCols));
 
  173 #if defined(KOKKOS_ENABLE_SERIAL) 
  174 template<
class ViewType,
 
  179 struct Fill<ViewType,
 
  186   fill (
const Kokkos::Serial& ,
 
  188         const ValueType& alpha,
 
  189         const IndexType numRows,
 
  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;
 
  197     typedef typename ViewType::non_const_value_type view_value_type;
 
  203     if (podType && X.span_is_contiguous () && alpha == ValueType (0.0)) {
 
  204       memsetWrapper (X.data (), 0, X.span () * 
sizeof (view_value_type));
 
  207       for (IndexType k = 0; k < numRows; ++k) {
 
  213 #endif // defined(KOKKOS_ENABLE_SERIAL) 
  215 #if defined(KOKKOS_ENABLE_SERIAL) 
  216 template<
class ViewType,
 
  221 struct Fill<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));
 
  267 #endif // defined(KOKKOS_ENABLE_SERIAL) 
  282 template<
class ViewType,
 
  285          class ExecutionSpace>
 
  287 fill (
const ExecutionSpace& execSpace,
 
  289       const ValueType& alpha,
 
  290       const IndexType numRows,
 
  291       const IndexType numCols)
 
  293   static_assert (std::is_integral<IndexType>::value,
 
  294                  "IndexType must be a built-in integer type.");
 
  295   typedef Impl::Fill<ViewType, ValueType, ExecutionSpace,
 
  296     IndexType, ViewType::Rank> impl_type;
 
  297   impl_type::fill (execSpace, X, alpha, numRows, numCols);
 
  300 template<
class ViewType,
 
  303          class ExecutionSpace>
 
  305 fill (
const ExecutionSpace& execSpace,
 
  307       const ValueType& alpha,
 
  308       const IndexType numRows,
 
  309       const IndexType numCols,
 
  310       const size_t whichVectors[])
 
  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.");
 
  317   for (IndexType k = 0; k < numCols; ++k) {
 
  318     const IndexType j = whichVectors[k];
 
  319     auto X_j = Kokkos::subview (X, Kokkos::ALL (), j);
 
  320     typedef decltype (X_j) one_d_view_type;
 
  321     typedef Impl::Fill<one_d_view_type, ValueType, ExecutionSpace,
 
  322       IndexType, 1> impl_type;
 
  323     impl_type::fill (execSpace, X_j, alpha, numRows, IndexType (1));
 
  331 #endif // TPETRA_DETAILS_FILL_HPP 
Type traits for Tpetra's BLAS wrappers; an implementation detail of Tpetra::MultiVector. 
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst. 
Implementation of ::Tpetra::Details::Blas::fill. 
Do BLAS libraries (all that are compliant with the BLAS Standard) support the given "scalar" (matrix ...