10#ifndef TPETRA_BLOCKVIEW_HPP
11#define TPETRA_BLOCKVIEW_HPP
21#include "TpetraCore_config.h"
22#include "KokkosKernels_ArithTraits.hpp"
23#include "Kokkos_Complex.hpp"
40template <
class ViewType1,
42 const int rank1 = ViewType1::rank>
59 static_assert(
static_cast<int>(ViewType1::rank) ==
static_cast<int>(ViewType2::rank),
60 "AbsMax: ViewType1 and ViewType2 must have the same rank.");
61 typedef typename std::remove_reference<
decltype(
Y(0, 0))>::type
STY;
62 static_assert(!std::is_const<STY>::value,
63 "AbsMax: The type of each entry of Y must be nonconst.");
64 typedef typename std::decay<
decltype(
X(0, 0))>::type
STX;
65 static_assert(std::is_same<STX, STY>::value,
66 "AbsMax: The type of each entry of X and Y must be the same.");
67 typedef KokkosKernels::ArithTraits<STY>
KAT;
69 const int numCols =
Y.extent(1);
71 for (
int j = 0;
j < numCols; ++
j) {
96 static_assert(
static_cast<int>(ViewType1::rank) ==
static_cast<int>(ViewType2::rank),
97 "AbsMax: ViewType1 and ViewType2 must have the same rank.");
99 typedef typename std::remove_reference<
decltype(
Y(0))>::type
STY;
100 static_assert(!std::is_const<STY>::value,
101 "AbsMax: The type of each entry of Y must be nonconst.");
102 typedef typename std::remove_const<
typename std::remove_reference<
decltype(
X(0))>::type>::type
STX;
103 static_assert(std::is_same<STX, STY>::value,
104 "AbsMax: The type of each entry of X and Y must be the same.");
105 typedef KokkosKernels::ArithTraits<STY>
KAT;
126template <
class ViewType1,
class ViewType2, const
int rank = ViewType1::rank>
129 static_assert(
static_cast<int>(ViewType1::rank) ==
static_cast<int>(ViewType2::rank),
130 "absMax: ViewType1 and ViewType2 must have the same rank.");
142 const int rank = ViewType::rank>
160#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
192 using x_value_type =
typename std::decay<
decltype(*
x.data())>::type;
195#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
207template <
class ViewType,
209 class IndexType = int,
210 const bool is_contiguous =
false,
211 const int rank = ViewType::rank>
254#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
266template <
class CoefficientType,
269 class IndexType = int,
270 const bool is_contiguous =
false,
271 const int rank = ViewType1::rank>
291 using KokkosKernels::ArithTraits;
292 static_assert(
static_cast<int>(ViewType1::rank) ==
static_cast<int>(ViewType2::rank),
293 "AXPY: x and y must have the same rank.");
316 using KokkosKernels::ArithTraits;
317 static_assert(
static_cast<int>(ViewType1::rank) ==
static_cast<int>(ViewType2::rank),
318 "AXPY: X and Y must have the same rank.");
341 using KokkosKernels::ArithTraits;
342 static_assert(
static_cast<int>(ViewType1::rank) ==
static_cast<int>(ViewType2::rank),
343 "AXPY: x and y must have the same rank.");
346 using x_value_type =
typename std::decay<
decltype(*
x.data())>::type;
347 using y_value_type =
typename std::decay<
decltype(*
y.data())>::type;
351#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
364template <
class ViewType1,
366 class IndexType = int,
367 const bool is_contiguous =
false,
368 const int rank = ViewType1::rank>
417 using x_value_type =
typename std::decay<
decltype(*
x.data())>::type;
418 using y_value_type =
typename std::decay<
decltype(*
y.data())>::type;
422#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
430template <
class VecType1,
434 class IndexType = int,
435 bool is_contiguous =
false,
436 class BlkLayoutType =
typename BlkType::array_layout>
438 static KOKKOS_INLINE_FUNCTION
void
439 run(
const CoeffType& alpha,
443 static_assert(VecType1::rank == 1,
"GEMV: VecType1 must have rank 1.");
444 static_assert(BlkType::rank == 2,
"GEMV: BlkType must have rank 2.");
445 static_assert(VecType2::rank == 1,
"GEMV: VecType2 must have rank 1.");
447 const IndexType numRows =
static_cast<IndexType
>(A.extent(0));
448 const IndexType numCols =
static_cast<IndexType
>(A.extent(1));
451 for (IndexType i = 0; i < numRows; ++i)
452 for (IndexType j = 0; j < numCols; ++j)
453 y(i) += alpha * A(i, j) * x(j);
457template <
class VecType1,
462struct GEMV<VecType1, BlkType, VecType2, CoeffType, IndexType, true, Kokkos::LayoutLeft> {
463 static KOKKOS_INLINE_FUNCTION
void
464 run(
const CoeffType& alpha,
468 static_assert(VecType1::rank == 1,
"GEMV: VecType1 must have rank 1.");
469 static_assert(BlkType::rank == 2,
"GEMV: BlkType must have rank 2.");
470 static_assert(VecType2::rank == 1,
"GEMV: VecType2 must have rank 1.");
472 using A_value_type =
typename std::decay<
decltype(A(0, 0))>::type;
473 using x_value_type =
typename std::decay<
decltype(x(0))>::type;
474 using y_value_type =
typename std::decay<
decltype(y(0))>::type;
476 const IndexType numRows =
static_cast<IndexType
>(A.extent(0));
477 const IndexType numCols =
static_cast<IndexType
>(A.extent(1));
479 const A_value_type* __restrict A_ptr(A.data());
480 const IndexType as1(A.stride(1));
481 const x_value_type* __restrict x_ptr(x.data());
482 y_value_type* __restrict y_ptr(y.data());
484#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
487 for (IndexType j = 0; j < numCols; ++j) {
488 const x_value_type x_at_j = alpha * x_ptr[j];
489 const A_value_type* __restrict A_at_j = A_ptr + j * as1;
490#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
493 for (IndexType i = 0; i < numRows; ++i)
494 y_ptr[i] += A_at_j[i] * x_at_j;
499template <
class VecType1,
504struct GEMV<VecType1, BlkType, VecType2, CoeffType, IndexType, true, Kokkos::LayoutRight> {
505 static KOKKOS_INLINE_FUNCTION
void
506 run(
const CoeffType& alpha,
510 static_assert(VecType1::rank == 1,
"GEMV: VecType1 must have rank 1.");
511 static_assert(BlkType::rank == 2,
"GEMV: BlkType must have rank 2.");
512 static_assert(VecType2::rank == 1,
"GEMV: VecType2 must have rank 1.");
514 using A_value_type =
typename std::decay<
decltype(A(0, 0))>::type;
515 using x_value_type =
typename std::decay<
decltype(x(0))>::type;
516 using y_value_type =
typename std::decay<
decltype(y(0))>::type;
518 const IndexType numRows =
static_cast<IndexType
>(A.extent(0));
519 const IndexType numCols =
static_cast<IndexType
>(A.extent(1));
521 const A_value_type* __restrict A_ptr(A.data());
522 const IndexType as0(A.stride(0));
523 const x_value_type* __restrict x_ptr(x.data());
524 y_value_type* __restrict y_ptr(y.data());
526#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
529 for (IndexType i = 0; i < numRows; ++i) {
530 y_value_type y_at_i(0);
531 const auto A_at_i = A_ptr + i * as0;
532#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
535 for (IndexType j = 0; j < numCols; ++j)
536 y_at_i += A_at_i[j] * x_ptr[j];
537 y_ptr[i] += alpha * y_at_i;
546template <
class ViewType,
547 class CoefficientType,
548 class IndexType = int,
549 const int rank = ViewType::rank>
550KOKKOS_INLINE_FUNCTION
void
552 using LayoutType =
typename ViewType::array_layout;
553 constexpr bool is_contiguous = (std::is_same<LayoutType, Kokkos::LayoutLeft>::value ||
554 std::is_same<LayoutType, Kokkos::LayoutRight>::value);
559template <
class ViewType,
561 class IndexType = int,
562 const int rank = ViewType::rank>
563KOKKOS_INLINE_FUNCTION
void
565 using LayoutType =
typename ViewType::array_layout;
566 constexpr bool is_contiguous = (std::is_same<LayoutType, Kokkos::LayoutLeft>::value ||
567 std::is_same<LayoutType, Kokkos::LayoutRight>::value);
576template <
class CoefficientType,
579 class IndexType = int,
580 const int rank = ViewType1::rank>
581KOKKOS_INLINE_FUNCTION
void
585 static_assert(
static_cast<int>(ViewType1::rank) ==
586 static_cast<int>(ViewType2::rank),
587 "AXPY: x and y must have the same rank.");
588 using LayoutType1 =
typename ViewType1::array_layout;
589 using LayoutType2 =
typename ViewType2::array_layout;
590 constexpr bool is_layout_same = std::is_same<LayoutType1, LayoutType2>::value;
591 constexpr bool is_x_contiguous = (std::is_same<LayoutType1, Kokkos::LayoutLeft>::value ||
592 std::is_same<LayoutType1, Kokkos::LayoutRight>::value);
605template <
class ViewType1,
607 class IndexType = int,
608 const int rank = ViewType1::rank>
609KOKKOS_INLINE_FUNCTION
void
611 static_assert(
static_cast<int>(ViewType1::rank) ==
612 static_cast<int>(ViewType2::rank),
613 "COPY: x and y must have the same rank.");
614 using LayoutType1 =
typename ViewType1::array_layout;
615 using LayoutType2 =
typename ViewType2::array_layout;
616 constexpr bool is_layout_same = std::is_same<LayoutType1, LayoutType2>::value;
617 constexpr bool is_x_contiguous = (std::is_same<LayoutType1, Kokkos::LayoutLeft>::value ||
618 std::is_same<LayoutType1, Kokkos::LayoutRight>::value);
634template <
class VecType1,
638 class IndexType =
int>
639KOKKOS_INLINE_FUNCTION
void
644 constexpr bool is_A_contiguous = (std::is_same<typename BlkType::array_layout, Kokkos::LayoutLeft>::value ||
645 std::is_same<typename BlkType::array_layout, Kokkos::LayoutRight>::value);
646 constexpr bool is_x_contiguous = (std::is_same<typename VecType1::array_layout, Kokkos::LayoutLeft>::value ||
647 std::is_same<typename VecType1::array_layout, Kokkos::LayoutRight>::value);
648 constexpr bool is_y_contiguous = (std::is_same<typename VecType2::array_layout, Kokkos::LayoutLeft>::value ||
649 std::is_same<typename VecType2::array_layout, Kokkos::LayoutRight>::value);
651 Impl::GEMV<VecType1, BlkType, VecType2, CoeffType, IndexType, is_contiguous>::run(
alpha,
A,
x,
y);
661template <
class ViewType1,
664 class CoefficientType,
665 class IndexType =
int>
666KOKKOS_INLINE_FUNCTION
void
675 static_assert(ViewType1::rank == 2,
"GEMM: A must have rank 2 (be a matrix).");
676 static_assert(ViewType2::rank == 2,
"GEMM: B must have rank 2 (be a matrix).");
677 static_assert(ViewType3::rank == 2,
"GEMM: C must have rank 2 (be a matrix).");
679 typedef typename std::remove_reference<
decltype(
A(0, 0))>::type
Scalar;
680 typedef KokkosKernels::ArithTraits<Scalar> STS;
793template <
class LittleBlockType,
794 class LittleVectorType>
795KOKKOS_INLINE_FUNCTION
void
798 typedef typename std::decay<
decltype(
ipiv(0))>::type
IndexType;
799 static_assert(std::is_integral<IndexType>::value,
800 "GETF2: The type of each entry of ipiv must be an integer type.");
801 typedef typename std::remove_reference<
decltype(
A(0, 0))>::type
Scalar;
802 static_assert(!std::is_const<Scalar>::value,
803 "GETF2: A must not be a const View (or LittleBlock).");
804 static_assert(!std::is_const<std::remove_reference<
decltype(
ipiv(0))>>::value,
805 "GETF2: ipiv must not be a const View (or LittleBlock).");
806 static_assert(LittleBlockType::rank == 2,
"GETF2: A must have rank 2 (be a matrix).");
807 typedef KokkosKernels::ArithTraits<Scalar> STS;
828 if (STS::abs(
A(
i,
j)) > STS::abs(
A(
jp,
j))) {
848 }
else if (
info == 0) {
866template <
class LittleBlockType,
867 class LittleIntVectorType,
868 class LittleScalarVectorType,
869 const int rank = LittleScalarVectorType::rank>
872 run(
const char mode[],
885 run(
const char mode[],
891 typedef typename std::decay<
decltype(
ipiv(0))>::type
IndexType;
894 static_assert(std::is_integral<IndexType>::value &&
895 std::is_signed<IndexType>::value,
896 "GETRS: The type of each entry of ipiv must be a signed integer.");
897 typedef typename std::decay<
decltype(
A(0, 0))>::type
Scalar;
898 static_assert(!std::is_const<std::remove_reference<
decltype(
B(0))>>::value,
899 "GETRS: B must not be a const View (or LittleBlock).");
900 static_assert(LittleBlockType::rank == 2,
"GETRS: A must have rank 2 (be a matrix).");
901 static_assert(LittleIntVectorType::rank == 1,
"GETRS: ipiv must have rank 1.");
902 static_assert(LittleScalarVectorType::rank == 1,
"GETRS: For this specialization, B must have rank 1.");
904 typedef KokkosKernels::ArithTraits<Scalar> STS;
925 if (
mode[0] ==
'n' ||
mode[0] ==
'N') {
957 else if (
mode[0] ==
't' ||
mode[0] ==
'T') {
962 else if (
mode[0] ==
'c' ||
mode[0] ==
'C') {
978 run(
const char mode[],
984 typedef typename std::decay<
decltype(
ipiv(0))>::type
IndexType;
985 static_assert(std::is_integral<IndexType>::value,
986 "GETRS: The type of each entry of ipiv must be an integer type.");
987 static_assert(!std::is_const<std::remove_reference<
decltype(
B(0))>>::value,
988 "GETRS: B must not be a const View (or LittleBlock).");
989 static_assert(LittleBlockType::rank == 2,
"GETRS: A must have rank 2 (be a matrix).");
990 static_assert(LittleIntVectorType::rank == 1,
"GETRS: ipiv must have rank 1.");
991 static_assert(LittleScalarVectorType::rank == 2,
"GETRS: For this specialization, B must have rank 2.");
1000 auto B_cur = Kokkos::subview(
B, Kokkos::ALL(),
rhs);
1041template <
class LittleBlockType,
1042 class LittleIntVectorType,
1043 class LittleScalarVectorType>
1044KOKKOS_INLINE_FUNCTION
void
1050 typedef typename std::decay<
decltype(
ipiv(0))>::type
IndexType;
1053 static_assert(std::is_integral<IndexType>::value &&
1054 std::is_signed<IndexType>::value,
1055 "GETRI: The type of each entry of ipiv must be a signed integer.");
1056 typedef typename std::remove_reference<
decltype(
A(0, 0))>::type
Scalar;
1057 static_assert(!std::is_const<std::remove_reference<
decltype(
A(0, 0))>>::value,
1058 "GETRI: A must not be a const View (or LittleBlock).");
1059 static_assert(!std::is_const<std::remove_reference<
decltype(
work(0))>>::value,
1060 "GETRI: work must not be a const View (or LittleBlock).");
1061 static_assert(LittleBlockType::rank == 2,
"GETRI: A must have rank 2 (be a matrix).");
1062 typedef KokkosKernels::ArithTraits<Scalar> STS;
1144template<
class LittleBlockType,
1145 class LittleVectorTypeX,
1146 class LittleVectorTypeY,
1147 class CoefficientType,
1148 class IndexType =
int>
1149KOKKOS_INLINE_FUNCTION
void
1150GEMV (
const char trans,
1151 const CoefficientType& alpha,
1152 const LittleBlockType& A,
1153 const LittleVectorTypeX& x,
1154 const CoefficientType& beta,
1155 const LittleVectorTypeY& y)
1161 typedef typename std::remove_reference<
decltype (y(0)) >::type y_value_type;
1162 const IndexType numRows =
static_cast<IndexType
> (A.extent (0));
1163 const IndexType numCols =
static_cast<IndexType
> (A.extent (1));
1167 for (IndexType i = 0; i < numRows; ++i) {
1172 for (IndexType i = 0; i < numRows; ++i) {
1173 y_value_type y_i = 0.0;
1174 for (IndexType j = 0; j < numCols; ++j) {
1175 y_i += A(i,j) * x(j);
1184 for (IndexType i = 0; i < numRows; ++i) {
1189 for (IndexType i = 0; i < numRows; ++i) {
1195 for (IndexType i = 0; i < numRows; ++i) {
1196 y_value_type y_i = beta * y(i);
1197 for (IndexType j = 0; j < numCols; ++j) {
1198 y_i += alpha * A(i,j) * x(j);
Struct that holds views of the contents of a CrsMatrix.
KOKKOS_INLINE_FUNCTION void absMax(const ViewType2 &Y, const ViewType1 &X)
Implementation of Tpetra's ABSMAX CombineMode for the small dense blocks in BlockCrsMatrix,...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
KOKKOS_INLINE_FUNCTION void GEMV(const CoeffType &alpha, const BlkType &A, const VecType1 &x, const VecType2 &y)
y := y + alpha * A * x (dense matrix-vector multiply)
KOKKOS_INLINE_FUNCTION void GETF2(const LittleBlockType &A, const LittleVectorType &ipiv, int &info)
Computes A = P*L*U.
KOKKOS_INLINE_FUNCTION void GETRI(const LittleBlockType &A, const LittleIntVectorType &ipiv, const LittleScalarVectorType &work, int &info)
Compute inverse of A, using result of GETRF or GETF2.
KOKKOS_INLINE_FUNCTION void AXPY(const CoefficientType &alpha, const ViewType1 &x, const ViewType2 &y)
y := y + alpha * x (dense vector or matrix update)
KOKKOS_INLINE_FUNCTION void FILL(const ViewType &x, const InputType &val)
Set every entry of x to val.
KOKKOS_INLINE_FUNCTION void SCAL(const CoefficientType &alpha, const ViewType &x)
x := alpha*x, where x is either rank 1 (a vector) or rank 2 (a matrix).
KOKKOS_INLINE_FUNCTION void COPY(const ViewType1 &x, const ViewType2 &y)
Deep copy x into y, where x and y are either rank 1 (vectors) or rank 2 (matrices) with the same dime...
KOKKOS_INLINE_FUNCTION void GETRS(const char mode[], const LittleBlockType &A, const LittleIntVectorType &ipiv, const LittleScalarVectorType &B, int &info)
Solve the linear system(s) AX=B, using the result of GETRF or GETF2.
KOKKOS_INLINE_FUNCTION void GEMM(const char transA[], const char transB[], const CoefficientType &alpha, const ViewType1 &A, const ViewType2 &B, const CoefficientType &beta, const ViewType3 &C)
Small dense matrix-matrix multiply: C := alpha*A*B + beta*C
@ ZERO
Replace old values with zero.
static KOKKOS_INLINE_FUNCTION void run(const CoefficientType &alpha, const ViewType1 &x, const ViewType2 &y)
y := y + alpha*x (rank-1 x and y, i.e., vectors)
static KOKKOS_INLINE_FUNCTION void run(const CoefficientType &alpha, const ViewType1 &X, const ViewType2 &Y)
Y := Y + alpha*X (rank-2 X and Y, i.e., matrices)
Implementation of Tpetra::AXPY function.
static KOKKOS_INLINE_FUNCTION void run(const ViewType2 &Y, const ViewType1 &X)
(*this)(i) := max(abs((*this)(i)), abs(X(i))) for all i.
static KOKKOS_INLINE_FUNCTION void run(const ViewType2 &Y, const ViewType1 &X)
(*this)(i,j) := max(abs((*this)(i,j)), abs(X(i,j))) for all (i,j).
Implementation of Tpetra's ABSMAX CombineMode for the small dense blocks in BlockCrsMatrix,...
static KOKKOS_INLINE_FUNCTION void run(const ViewType1 &x, const ViewType2 &y)
y := x (rank-1 x and y, i.e., vectors)
static KOKKOS_INLINE_FUNCTION void run(const ViewType1 &X, const ViewType2 &Y)
Y := X (rank-2 X and Y, i.e., matrices)
Implementation of Tpetra::COPY function.
Implementation of Tpetra::FILL function.
Computes the solution to Ax=b.
static KOKKOS_INLINE_FUNCTION void run(const CoefficientType &alpha, const ViewType &x)
x := alpha*x (rank-1 x, i.e., a vector)
static KOKKOS_INLINE_FUNCTION void run(const CoefficientType &alpha, const ViewType &A)
A := alpha*A (rank-2 A, i.e., a matrix)
Implementation of Tpetra::SCAL function.