46 #include "Stokhos_Ensemble_Sizes.hpp"
51 #include "Kokkos_Core.hpp"
54 template<
typename IntType >
61 return k + N * ( j + N * i );
64 template <
typename ordinal >
67 std::vector< std::vector<ordinal> >& graph )
69 graph.resize( N * N * N, std::vector<ordinal>() );
73 for (
int i = 0; i < (
int) N; ++i ) {
74 for (
int j = 0;
j < (
int) N; ++
j ) {
75 for (
int k = 0; k < (
int) N; ++k ) {
79 graph[row].reserve(27);
81 for (
int ii = -1; ii < 2; ++ii ) {
82 for (
int jj = -1; jj < 2; ++jj ) {
83 for (
int kk = -1; kk < 2; ++kk ) {
84 if ( 0 <= i + ii && i + ii < (
int) N &&
85 0 <=
j + jj &&
j + jj < (
int) N &&
86 0 <= k + kk && k + kk < (
int) N ) {
89 graph[row].push_back(col);
92 total += graph[row].size();
98 template <
typename scalar,
typename ordinal>
101 const ordinal nStoch,
102 const ordinal iRowFEM,
103 const ordinal iColFEM,
104 const ordinal iStoch )
106 const scalar A_fem = ( 10.0 + scalar(iRowFEM) / scalar(nFEM) ) +
107 ( 5.0 + scalar(iColFEM) / scalar(nFEM) );
109 const scalar A_stoch = ( 1.0 + scalar(iStoch) / scalar(nStoch) );
111 return A_fem + A_stoch;
115 template <
typename scalar,
typename ordinal>
118 const ordinal nStoch,
119 const ordinal iColFEM,
120 const ordinal iStoch )
122 const scalar X_fem = 100.0 + scalar(iColFEM) / scalar(nFEM);
123 const scalar X_stoch = 1.0 + scalar(iStoch) / scalar(nStoch);
124 return X_fem + X_stoch;
130 template <>
struct ScalarTol<float> {
static float tol() {
return 1e-4; } };
134 template <
typename array_type,
typename scalar_type>
136 const array_type& y_exp,
141 typedef typename array_type::size_type size_type;
147 size_type num_rows = y.extent(0);
148 size_type num_cols = y.extent(1);
150 for (size_type i=0; i<num_rows; ++i) {
151 for (size_type
j=0;
j<num_cols; ++
j) {
155 out <<
"y_expected(" << i <<
"," <<
j <<
") - "
156 <<
"y(" << i <<
"," <<
j <<
") = " << hy_exp(i,
j)
157 <<
" - " << hy(i,
j) <<
" == "
158 << diff <<
" < " << tol <<
" : ";
164 success = success && s;
172 template <
typename MatrixType>
177 typedef typename MatrixType::StaticCrsGraphType matrix_graph_type;
178 typedef typename MatrixType::values_type matrix_values_type;
180 std::vector< std::vector<ordinal_type> > graph(nrow);
181 for (ordinal_type i=0; i<nrow; ++i)
182 graph[i] = std::vector<ordinal_type>(1, i);
183 ordinal_type graph_length = nrow;
185 matrix_graph_type matrix_graph =
186 Kokkos::create_staticcrsgraph<matrix_graph_type>(
"graph", graph);
187 matrix_values_type matrix_values =
188 matrix_values_type(
"values", graph_length, mp_vector_size);
190 MatrixType matrix(
"matrix", nrow, matrix_values, matrix_graph);
199 template <
typename MatrixType>
210 KOKKOS_INLINE_FUNCTION
215 m_matrix.replaceValues(row, &col, 1, &val,
false,
true);
219 static void apply(
const MatrixType matrix) {
225 static bool check(
const MatrixType matrix,
227 typedef typename MatrixType::values_type matrix_values_type;
228 typename matrix_values_type::HostMirror host_matrix_values =
236 "matrix_values(row)",
240 success = success && s;
247 template <
typename MatrixType>
258 KOKKOS_INLINE_FUNCTION
263 m_matrix.sumIntoValues(row, &col, 1, &val,
false,
true);
267 static void apply(
const MatrixType matrix) {
273 static bool check(
const MatrixType matrix,
275 typedef typename MatrixType::values_type matrix_values_type;
276 typename matrix_values_type::HostMirror host_matrix_values =
284 "matrix_values(row)",
288 success = success && s;
296 template <
typename MatrixType>
307 KOKKOS_INLINE_FUNCTION
312 m_matrix.sumIntoValues(row, &col, 1, &val,
false,
true);
316 static void apply(
const MatrixType matrix) {
322 static bool check(
const MatrixType matrix,
324 typedef typename MatrixType::values_type matrix_values_type;
325 typename matrix_values_type::HostMirror host_matrix_values =
334 val =
value_type( vec_size, nrow*(nrow-1)/2 );
338 "matrix_values(row)",
342 success = success && s;
351 Kokkos_CrsMatrix_MP, ReplaceValues, MatrixScalar )
355 typedef KokkosSparse::CrsMatrix<MatrixScalar,Ordinal,Device> Matrix;
359 Matrix matrix = buildDiagonalMatrix<Matrix>(nrow,
VectorSize);
363 kernel::apply(matrix);
366 success = kernel::check(matrix, out);
370 Kokkos_CrsMatrix_MP, SumIntoValues, MatrixScalar )
374 typedef KokkosSparse::CrsMatrix<MatrixScalar,Ordinal,Device> Matrix;
378 Matrix matrix = buildDiagonalMatrix<Matrix>(nrow,
VectorSize);
382 kernel::apply(matrix);
385 success = kernel::check(matrix, out);
389 Kokkos_CrsMatrix_MP, SumIntoValuesAtomic, MatrixScalar )
393 typedef KokkosSparse::CrsMatrix<MatrixScalar,Ordinal,Device> Matrix;
397 Matrix matrix = buildDiagonalMatrix<Matrix>(nrow,
VectorSize);
401 kernel::apply(matrix);
404 success = kernel::check(matrix, out);
411 template <
class ViewType,
class OrdinalType,
size_t I >
413 static ViewType
create_view(
const std::string & name,
const OrdinalType & fem_length,
const OrdinalType & stoch_length ) {
414 return ViewType(name, fem_length, stoch_length);
418 template <
class ViewType,
class OrdinalType>
420 static ViewType
create_view(
const std::string & name,
const OrdinalType & fem_length,
const OrdinalType & stoch_length ) {
422 return ViewType(name, fem_length);
426 template <
class ViewType,
class OrdinalType>
428 static ViewType
create_view(
const std::string & name,
const OrdinalType & fem_length,
const OrdinalType & stoch_length ) {
431 return ViewType(name);
436 template <
typename VectorType,
typename Multiply>
439 KokkosSparse::DeviceConfig dev_config,
440 Multiply multiply_op,
447 typedef Kokkos::LayoutRight Layout;
448 typedef Kokkos::View< VectorType*, Layout, execution_space > block_vector_type;
449 typedef KokkosSparse::CrsMatrix< VectorType, ordinal_type, execution_space > block_matrix_type;
450 typedef typename block_matrix_type::StaticCrsGraphType matrix_graph_type;
451 typedef typename block_matrix_type::values_type matrix_values_type;
455 storage_type::is_static && storage_type::static_size != stoch_length,
457 "Static storage size must equal ensemble size");
460 ordinal_type fem_length = nGrid * nGrid * nGrid;
461 std::vector< std::vector<ordinal_type> > fem_graph;
474 block_vector_type x =
475 block_vector_type(
"x", fem_length, stoch_length);
476 block_vector_type y =
477 block_vector_type(
"y", fem_length, stoch_length);
483 typename block_vector_type::HostMirror::array_type hax = hx ;
484 typename block_vector_type::HostMirror::array_type hay = hy ;
486 for (ordinal_type iRowFEM=0; iRowFEM<fem_length; ++iRowFEM) {
487 for (ordinal_type iRowStoch=0; iRowStoch<stoch_length; ++iRowStoch) {
488 hax(iRowFEM,iRowStoch) =
489 generate_vector_coefficient<scalar_type>(
490 fem_length, stoch_length, iRowFEM, iRowStoch );
491 hay(iRowFEM,iRowStoch) = 0.0;
501 matrix_graph_type matrix_graph =
502 Kokkos::create_staticcrsgraph<matrix_graph_type>(
503 std::string(
"test crs graph"), fem_graph);
508 matrix_values_type matrix_values =
509 matrix_values_type(
"matrix", fem_graph_length, stoch_length);
510 block_matrix_type matrix(
511 "block_matrix", fem_length, matrix_values, matrix_graph);
512 matrix.dev_config = dev_config;
514 typename matrix_values_type::HostMirror hM =
517 typename matrix_values_type::HostMirror::array_type haM = hM ;
519 for (ordinal_type iRowFEM=0, iEntryFEM=0; iRowFEM<fem_length; ++iRowFEM) {
520 const ordinal_type row_size = fem_graph[iRowFEM].size();
521 for (ordinal_type iRowEntryFEM=0; iRowEntryFEM<row_size;
522 ++iRowEntryFEM, ++iEntryFEM) {
523 const ordinal_type iColFEM = fem_graph[iRowFEM][iRowEntryFEM];
525 for (ordinal_type k=0; k<stoch_length; ++k) {
527 generate_matrix_coefficient<scalar_type>(
528 fem_length, stoch_length, iRowFEM, iColFEM, k);
538 multiply_op( matrix, x, y );
543 typedef typename block_vector_type::array_type array_type;
544 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
545 array_type ay_expected =
546 array_type(
"ay_expected", fem_length, stoch_length);
548 array_type ay_expected =
551 typename array_type::HostMirror hay_expected =
553 for (ordinal_type iRowFEM=0, iEntryFEM=0; iRowFEM<fem_length; ++iRowFEM) {
554 const ordinal_type row_size = fem_graph[iRowFEM].size();
555 for (ordinal_type iRowEntryFEM=0; iRowEntryFEM<row_size;
556 ++iRowEntryFEM, ++iEntryFEM) {
557 const ordinal_type iColFEM = fem_graph[iRowFEM][iRowEntryFEM];
558 for (ordinal_type k=0; k<stoch_length; ++k) {
559 hay_expected(iRowFEM, k) +=
560 generate_matrix_coefficient<scalar_type>(
561 fem_length, stoch_length, iRowFEM, iColFEM, k) *
562 generate_vector_coefficient<scalar_type>(
563 fem_length, stoch_length, iColFEM, k );
572 typename block_vector_type::array_type ay = y;
581 template <
typename Matrix,
typename InputVector,
typename OutputVector>
583 const InputVector& x,
584 OutputVector& y)
const {
589 template <
typename Tag>
594 template <
typename Matrix,
typename InputVector,
typename OutputVector>
596 const InputVector& x,
597 OutputVector& y)
const {
605 #define CRSMATRIX_MP_VECTOR_TESTS_MATRIXSCALAR( SCALAR ) \
606 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( \
607 Kokkos_CrsMatrix_MP, ReplaceValues, SCALAR ) \
608 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( \
609 Kokkos_CrsMatrix_MP, SumIntoValues, SCALAR ) \
610 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( \
611 Kokkos_CrsMatrix_MP, SumIntoValuesAtomic, SCALAR )
613 #define CRSMATRIX_MP_VECTOR_TESTS_STORAGE( STORAGE ) \
614 typedef Sacado::MP::Vector<STORAGE> MP_Vector_ ## STORAGE; \
615 CRSMATRIX_MP_VECTOR_TESTS_MATRIXSCALAR( MP_Vector_ ## STORAGE )
617 #define CRSMATRIX_MP_VECTOR_TESTS_ORDINAL_SCALAR_DEVICE( ORDINAL, SCALAR, DEVICE ) \
618 typedef Stokhos::StaticFixedStorage<ORDINAL,SCALAR,VectorSize,DEVICE> SFS; \
619 typedef Stokhos::DynamicStorage<ORDINAL,SCALAR,DEVICE> DS; \
620 CRSMATRIX_MP_VECTOR_TESTS_STORAGE( SFS ) \
621 CRSMATRIX_MP_VECTOR_TESTS_STORAGE( DS )
623 #define CRSMATRIX_MP_VECTOR_TESTS_DEVICE( DEVICE ) \
624 CRSMATRIX_MP_VECTOR_TESTS_ORDINAL_SCALAR_DEVICE( int, double, DEVICE )
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
MatrixType::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
Stokhos::StandardStorage< int, double > storage_type
MatrixType::size_type size_type
const MatrixType m_matrix
MatrixType::value_type value_type
bool compareVecs(const VectorType1 &a1, const std::string &a1_name, const VectorType2 &a2, const std::string &a2_name, const ValueType &rel_tol, const ValueType &abs_tol, Teuchos::FancyOStream &out)
ordinal generate_fem_graph(ordinal N, std::vector< std::vector< ordinal > > &graph)
static bool check(const MatrixType matrix, Teuchos::FancyOStream &out)
Stokhos_MV_Multiply_Op< Stokhos::DefaultMultiply > DefaultMultiply
Kokkos::DefaultExecutionSpace execution_space
MatrixType::execution_space execution_space
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
MatrixType buildDiagonalMatrix(typename MatrixType::ordinal_type nrow, typename MatrixType::ordinal_type mp_vector_size)
MatrixType::ordinal_type ordinal_type
static void apply(const MatrixType matrix)
scalar generate_vector_coefficient(const ordinal nFEM, const ordinal nStoch, const ordinal iColFEM, const ordinal iStoch)
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< is_view_uq_pce< View< T, P...> >::value, unsigned >::type dimension_scalar(const View< T, P...> &view)
IntType map_fem_graph_coord(const IntType &N, const IntType &i, const IntType &j, const IntType &k)
MatrixType::value_type value_type
MatrixType::ordinal_type ordinal_type
ReplaceDiagonalValuesKernel(const MatrixType matrix)
MatrixType::size_type size_type
void multiply(const CrsMatrix< MatrixValue, Device, Layout > &A, const InputMultiVectorType &x, OutputMultiVectorType &y, const std::vector< OrdinalType > &col_indices, SingleColumnMultivectorMultiply)
static bool check(const MatrixType matrix, Teuchos::FancyOStream &out)
AddDiagonalValuesAtomicKernel(const MatrixType matrix)
Kokkos_MV_Multiply_Op KokkosMultiply
MatrixType::execution_space execution_space
void operator()(const Matrix &A, const InputVector &x, OutputVector &y) const
static void apply(const MatrixType matrix)
MatrixType::size_type size_type
static bool check(const MatrixType matrix, Teuchos::FancyOStream &out)
Stokhos_MV_Multiply_Op(const Tag &tg=Tag())
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
static ViewType create_view(const std::string &name, const OrdinalType &fem_length, const OrdinalType &stoch_length)
static ViewType create_view(const std::string &name, const OrdinalType &fem_length, const OrdinalType &stoch_length)
static ViewType create_view(const std::string &name, const OrdinalType &fem_length, const OrdinalType &stoch_length)
scalar generate_matrix_coefficient(const ordinal nFEM, const ordinal nStoch, const ordinal iRowFEM, const ordinal iColFEM, const ordinal iStoch)
const MatrixType m_matrix
KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL(Kokkos_CrsMatrix_MP, ReplaceValues, MatrixScalar)
AddDiagonalValuesKernel(const MatrixType matrix)
MatrixType::ordinal_type ordinal_type
const MatrixType m_matrix
void operator()(const Matrix &A, const InputVector &x, OutputVector &y) const
bool compare_rank_2_views(const array_type &y, const array_type &y_exp, const scalar_type rel_tol, const scalar_type abs_tol, Teuchos::FancyOStream &out)
static void apply(const MatrixType matrix)
const unsigned VectorSize
bool test_embedded_vector(const typename VectorType::ordinal_type nGrid, const typename VectorType::ordinal_type stoch_length, KokkosSparse::DeviceConfig dev_config, Multiply multiply_op, Teuchos::FancyOStream &out)
MatrixType::value_type value_type
std::enable_if< Kokkos::is_view_uq_pce< Kokkos::View< InputType, InputP... > >::value &&Kokkos::is_view_uq_pce< Kokkos::View< OutputType, OutputP... > >::value >::type spmv(const char mode[], const AlphaType &a, const MatrixType &A, const Kokkos::View< InputType, InputP... > &x, const BetaType &b, const Kokkos::View< OutputType, OutputP... > &y, const RANK_ONE)
Stokhos::CrsMatrix< ValueType, Device, Layout >::HostMirror create_mirror_view(const Stokhos::CrsMatrix< ValueType, Device, Layout > &A)