44#ifndef TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_DIST_OBJECT_KERNELS_UQ_PCE_HPP
45#define TPETRA_KOKKOS_REFACTOR_DETAILS_MULTI_VECTOR_DIST_OBJECT_KERNELS_UQ_PCE_HPP
51#include "Tpetra_KokkosRefactor_Details_MultiVectorDistObjectKernels.hpp"
62 template <
typename DstView,
typename SrcView,
typename IdxView>
64 DstView, SrcView, IdxView,
65 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
66 Kokkos::is_view_uq_pce<SrcView>::value >::type >
69 typedef typename execution_space::size_type
size_type;
84 KOKKOS_INLINE_FUNCTION
89 KOKKOS_INLINE_FUNCTION
92 dst(k).fastAccessCoeff(i) =
src(
idx(k),
col).fastAccessCoeff(i);
100 Kokkos::parallel_for(
106 template <
typename DstView,
typename SrcView,
typename IdxView>
108 DstView, SrcView, IdxView,
109 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
110 Kokkos::is_view_uq_pce<SrcView>::value >::type >
128 KOKKOS_INLINE_FUNCTION
130 const typename IdxView::value_type localRow =
idx(k);
131 const size_t offset = k*
numCols;
136 KOKKOS_INLINE_FUNCTION
138 const typename IdxView::value_type localRow =
idx(k);
139 const size_t offset = k*
numCols;
142 dst(offset +
j).fastAccessCoeff(i) =
143 src(localRow,
j).fastAccessCoeff(i);
151 Kokkos::parallel_for(
157 template <
typename DstView,
typename SrcView,
typename IdxView,
160 DstView, SrcView, IdxView, ColView,
161 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
162 Kokkos::is_view_uq_pce<SrcView>::value >::type >
182 KOKKOS_INLINE_FUNCTION
184 const typename IdxView::value_type localRow =
idx(k);
185 const size_t offset = k*
numCols;
190 KOKKOS_INLINE_FUNCTION
192 const typename IdxView::value_type localRow =
idx(k);
193 const size_t offset = k*
numCols;
196 dst(offset +
j).fastAccessCoeff(i) =
197 src(localRow,
col(
j)).fastAccessCoeff(i);
206 Kokkos::parallel_for(
212 template <
typename ExecutionSpace,
218 ExecutionSpace, DstView, SrcView, IdxView, Op,
219 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
220 Kokkos::is_view_uq_pce<SrcView>::value >::type >
238 const size_t numCols_) :
246 template <
class TagType>
247 KOKKOS_INLINE_FUNCTION
void
250 const typename IdxView::value_type localRow =
idx(k);
251 const size_t offset = k*
numCols;
257 template <
class TagType>
258 KOKKOS_INLINE_FUNCTION
void
261 const typename IdxView::value_type localRow =
idx(k);
262 const size_t offset = k*
numCols;
278 const bool use_atomic_updates)
280 if (use_atomic_updates) {
282 (
"Tpetra::MultiVector (Sacado::UQ::PCE) unpack (constant stride)",
289 (
"Tpetra::MultiVector (Sacado::UQ::PCE) unpack (constant stride)",
297 template <
typename ExecutionSpace,
304 ExecutionSpace, DstView, SrcView, IdxView, ColView, Op,
305 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
306 Kokkos::is_view_uq_pce<SrcView>::value >::type>
329 template <
class TagType>
330 KOKKOS_INLINE_FUNCTION
void
332 const typename IdxView::value_type localRow =
idx(k);
333 const size_t offset = k*
numCols;
338 template <
class TagType>
339 KOKKOS_INLINE_FUNCTION
void
341 const typename IdxView::value_type localRow =
idx(k);
342 const size_t offset = k*
numCols;
357 const bool use_atomic_updates)
359 if (use_atomic_updates) {
361 (
"Tpetra::MultiVector unpack (Sacado::UQ::PCE) (nonconstant stride)",
369 (
"Tpetra::MultiVector unpack (Sacado::UQ::PCE) (nonconstant stride)",
378 template <
typename DstView,
typename SrcView,
379 typename DstIdxView,
typename SrcIdxView,
typename Op>
381 DstView, SrcView, DstIdxView, SrcIdxView, Op,
382 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
383 Kokkos::is_view_uq_pce<SrcView>::value >::type>
399 const DstIdxView& dst_idx_,
400 const SrcIdxView& src_idx_,
401 size_t numCols_,
const Op& op_) :
405 KOKKOS_INLINE_FUNCTION
407 const typename DstIdxView::value_type toRow =
dst_idx(k);
408 const typename SrcIdxView::value_type fromRow =
src_idx(k);
414 KOKKOS_INLINE_FUNCTION
416 const typename DstIdxView::value_type toRow =
dst_idx(k);
417 const typename SrcIdxView::value_type fromRow =
src_idx(k);
425 template <
typename ExecutionSpace>
426 static void permute(
const ExecutionSpace &space,
434 Kokkos::parallel_for(
440 template <
typename DstView,
typename SrcView,
441 typename DstIdxView,
typename SrcIdxView,
442 typename DstColView,
typename SrcColView,
typename Op>
444 DstView, SrcView, DstIdxView, SrcIdxView, DstColView, SrcColView, Op,
445 typename
std::enable_if< Kokkos::is_view_uq_pce<DstView>::value &&
446 Kokkos::is_view_uq_pce<SrcView>::value >::type >
464 const DstIdxView& dst_idx_,
465 const SrcIdxView& src_idx_,
466 const DstColView& dst_col_,
467 const SrcColView& src_col_,
474 KOKKOS_INLINE_FUNCTION
476 const typename DstIdxView::value_type toRow =
dst_idx(k);
477 const typename SrcIdxView::value_type fromRow =
src_idx(k);
483 KOKKOS_INLINE_FUNCTION
485 const typename DstIdxView::value_type toRow =
dst_idx(k);
486 const typename SrcIdxView::value_type fromRow =
src_idx(k);
494 template <
typename ExecutionSpace>
495 static void permute(
const ExecutionSpace &space,
505 Kokkos::parallel_for(
expr1 expr1 expr1 expr2 expr1 expr1 c expr2 expr1 c fastAccessCoeff(j) - expr2.val(j)
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< is_view_uq_pce< View< T, P... > >::value, unsigned >::type dimension_scalar(const View< T, P... > &view)
Team-based parallel work configuration for Sacado::MP::Vector.
static const unsigned BlockSize
PackArrayMultiColumnVariableStride(const DstView &dst_, const SrcView &src_, const IdxView &idx_, const ColView &col_, size_t numCols_)
KOKKOS_INLINE_FUNCTION void operator()(const size_type k) const
DstView::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(const size_type k, const size_type tidx) const
static void pack(const DstView &dst, const SrcView &src, const IdxView &idx, const ColView &col, size_t numCols, const execution_space &space)
execution_space::size_type size_type
PackArrayMultiColumn(const DstView &dst_, const SrcView &src_, const IdxView &idx_, size_t numCols_)
KOKKOS_INLINE_FUNCTION void operator()(const size_type k, const size_type tidx) const
KOKKOS_INLINE_FUNCTION void operator()(const size_type k) const
static const unsigned BlockSize
DstView::execution_space execution_space
execution_space::size_type size_type
static void pack(const DstView &dst, const SrcView &src, const IdxView &idx, size_t numCols, const execution_space &space)
execution_space::size_type size_type
static const unsigned BlockSize
DstView::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(const size_type k) const
static void pack(const DstView &dst, const SrcView &src, const IdxView &idx, size_t col, const execution_space &space)
KOKKOS_INLINE_FUNCTION void operator()(const size_type k, const size_type tidx) const
PackArraySingleColumn(const DstView &dst_, const SrcView &src_, const IdxView &idx_, size_t col_)
execution_space::size_type size_type
KOKKOS_INLINE_FUNCTION void operator()(const size_type k, const size_type tidx) const
KOKKOS_INLINE_FUNCTION void operator()(const size_type k) const
DstView::execution_space execution_space
static const unsigned BlockSize
static void permute(const ExecutionSpace &space, const DstView &dst, const SrcView &src, const DstIdxView &dst_idx, const SrcIdxView &src_idx, const DstColView &dst_col, const SrcColView &src_col, size_t numCols, const Op &op)
PermuteArrayMultiColumnVariableStride(const DstView &dst_, const SrcView &src_, const DstIdxView &dst_idx_, const SrcIdxView &src_idx_, const DstColView &dst_col_, const SrcColView &src_col_, size_t numCols_, const Op &op_)
execution_space::size_type size_type
static const unsigned BlockSize
PermuteArrayMultiColumn(const DstView &dst_, const SrcView &src_, const DstIdxView &dst_idx_, const SrcIdxView &src_idx_, size_t numCols_, const Op &op_)
DstView::execution_space execution_space
static void permute(const ExecutionSpace &space, const DstView &dst, const SrcView &src, const DstIdxView &dst_idx, const SrcIdxView &src_idx, size_t numCols, const Op &op)
KOKKOS_INLINE_FUNCTION void operator()(const size_type k) const
KOKKOS_INLINE_FUNCTION void operator()(const size_type k, const size_type tidx) const
KOKKOS_INLINE_FUNCTION void operator()(TagType tag, const size_type k) const
ExecutionSpace::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(TagType tag, const size_type k, const size_type tidx) const
execution_space::size_type size_type
static void unpack(const ExecutionSpace &execSpace, const DstView &dst, const SrcView &src, const IdxView &idx, const ColView &col, const Op &op, const size_t numCols, const bool use_atomic_updates)
static const unsigned BlockSize
UnpackArrayMultiColumnVariableStride(const ExecutionSpace &, const DstView &dst_, const SrcView &src_, const IdxView &idx_, const ColView &col_, const Op &op_, size_t numCols_)
static const unsigned BlockSize
ExecutionSpace::execution_space execution_space
execution_space::size_type size_type
static void unpack(const ExecutionSpace &execSpace, const DstView &dst, const SrcView &src, const IdxView &idx, const Op &op, const size_t numCols, const bool use_atomic_updates)
UnpackArrayMultiColumn(const ExecutionSpace &, const DstView &dst_, const SrcView &src_, const IdxView &idx_, const Op &op_, const size_t numCols_)