You signed in with another tab or window. Reload to refresh your session.You signed out in another tab or window. Reload to refresh your session.You switched accounts on another tab or window. Reload to refresh your session.Dismiss alert
Kokkos-Kernels does not compile when enabling the Trilinos option TPL_ENABLE_MKL=ON.
There are two issues. First the file src/sparse/impl/KokkosSparse_spgemm_mkl_impl.hpp generates errors such as
/home/etphipp/Trilinos/Trilinos/packages/kokkos-kernels/src/sparse/impl/KokkosSparse_spgemm_mkl_impl.hpp:212:9: error:
unknown type name 'sparse_index_base_t'; did you mean
'KokkosBatched::sparse_index_base_t'?
sparse_index_base_t c_indexing;
^~~~~~~~~~~~~~~~~~~
KokkosBatched::sparse_index_base_t
/home/software/intel/compilers_and_libraries_2019.3.199/linux/mkl/include/mkl_spblas.h:658:7: note:
'KokkosBatched::sparse_index_base_t' declared here
} sparse_index_base_t;
/home/etphipp/Trilinos/Trilinos/packages/kokkos-kernels/src/sparse/impl/KokkosSparse_spgemm_mkl_impl.hpp:216:13: error:
use of undeclared identifier 'SPARSE_STATUS_SUCCESS'; did you mean
'KokkosBatched::SPARSE_STATUS_SUCCESS'?
if (SPARSE_STATUS_SUCCESS !=
^~~~~~~~~~~~~~~~~~~~~
KokkosBatched::SPARSE_STATUS_SUCCESS
/home/software/intel/compilers_and_libraries_2019.3.199/linux/mkl/include/mkl_spblas.h:624:9: note:
'KokkosBatched::SPARSE_STATUS_SUCCESS' declared here
SPARSE_STATUS_SUCCESS = 0, /* the operation was suc...
Looking through that file, I thought I was seeing these because I had KokkosKernels_ENABLE_TPL_MKL=OFF (which is the default, regardless of TPL_ENABLE_MKL), but even turning that on doesn't fix it.
Second, files such as src/batched/KokkosBatched_Gemm_Serial_Impl.hpp and src/batched/KokkosBatched_Trsm_Serial_Impl.hpp also include MKL-specific specializations that generate errors if deprecated code in Kokkos is turned off, e.g.,
/home/etphipp/Trilinos/Trilinos/packages/kokkos-kernels/src/batched/KokkosBatched_Gemm_Serial_Impl.hpp:49:13: error: cannot refer to type member 'dimension' in 'const Kokkos::View<KokkosBatched::Vector<KokkosBatched::SIMD<double>, 4> **, Kokkos::LayoutStride, Kokkos::Device<Kokkos::OpenMP, Kokkos::HostSpace>, Kokkos::MemoryTraits<1> >' with '.'
m = C.dimension(0),
Replacing those lines with extent(0) fixes those errors. However those specializations appear to be compiled even with KokkosKernels_ENABLE_TPL_MKL=OFF, which I'm guessing they shouldn't be.
The text was updated successfully, but these errors were encountered:
@trilinos/kokkos-kernels
Kokkos-Kernels does not compile when enabling the Trilinos option TPL_ENABLE_MKL=ON.
There are two issues. First the file src/sparse/impl/KokkosSparse_spgemm_mkl_impl.hpp generates errors such as
Looking through that file, I thought I was seeing these because I had KokkosKernels_ENABLE_TPL_MKL=OFF (which is the default, regardless of TPL_ENABLE_MKL), but even turning that on doesn't fix it.
Second, files such as src/batched/KokkosBatched_Gemm_Serial_Impl.hpp and src/batched/KokkosBatched_Trsm_Serial_Impl.hpp also include MKL-specific specializations that generate errors if deprecated code in Kokkos is turned off, e.g.,
Replacing those lines with
extent(0)
fixes those errors. However those specializations appear to be compiled even with KokkosKernels_ENABLE_TPL_MKL=OFF, which I'm guessing they shouldn't be.The text was updated successfully, but these errors were encountered: