[Trilinos-Users] [EXTERNAL] Trilinos-Kokkos with CUDA

Phipps, Eric T etphipp at sandia.gov
Tue Nov 1 11:29:53 EDT 2016


Are you using nvcc_wrapper as your compiler?

-Eric

From: Trilinos-Users <trilinos-users-bounces at trilinos.org> on behalf of Ashesh Chattopadhyay <ashesh6810 at gmail.com>
Date: Tuesday, November 1, 2016 at 9:25 AM
To: "Trilinos-Users at trilinos.org" <Trilinos-Users at trilinos.org>
Subject: [EXTERNAL] [Trilinos-Users] Trilinos-Kokkos with CUDA

Hi,

I've built Trilinos with CUDA support and am trying to compile a simple program with CUDA.

 #include <limits>
#include <cstdio>
#include <cstdlib>
#include <cstring>
#include <sys/time.h>

#include <Kokkos_Core.hpp>

void checkSizes(int &N, int &M, int &S, int &nrepeat);

int main(int argc, char* argv[])
{

  int N = 1024 ;       // number of rows 2^12
  int M = 1024 ;       // number of columns 2^10
  int S = 1024*1024 ;      // total size 2^22
  int nrepeat = 100 ;    // number of repeats of the test






  Kokkos::initialize(argc,argv);

  // typedef Kokkos::Serial   ExecSpace ;
  // typedef Kokkos::Threads  ExecSpace ;
  // typedef Kokkos::OpenMP   ExecSpace ;
  typedef Kokkos::Cuda        ExecSpace ;

  //EXERCISE: Choose device memory space
  // typedef Kokkos::HostSpace    MemSpace;
  // typedef Kokkos::OpenMP       MemSpace;
  typedef Kokkos::CudaSpace       MemSpace;
  // typedef Kokkos::CudaUVMSpace MemSpace;

  typedef Kokkos::LayoutLeft   Layout ;
  // typedef Kokkos::LayoutRight  Layout ;

  typedef Kokkos::RangePolicy<ExecSpace> range_policy ;

  // Allocate y, x vectors and Matrix A:
  // Device
  typedef Kokkos::View<double*, Layout, MemSpace>   ViewVectorType;
  typedef Kokkos::View<double**, Layout, MemSpace>   ViewMatrixType;
  ViewVectorType devy("devy", N);
  ViewVectorType devx("devx", M);
  ViewMatrixType devA("devA", N, M);

  //Host mirror
  ViewVectorType::HostMirror y =  Kokkos::create_mirror_view(devy);
  ViewVectorType::HostMirror x =  Kokkos::create_mirror_view(devx);
  ViewMatrixType::HostMirror A =  Kokkos::create_mirror_view(devA);

  // Initialize y vector on host
  for (int i = 0; i < N; ++i) {
    y( i ) = 1;
  }

  // Initialize x vector on host
  for (int i = 0; i < M; ++i) {
    x( i ) = 1;
  }

  // Initialize A matrix, note 2D indexing computation on host
  for (int j = 0; j < N; ++j) {
    for ( int i = 0 ; i < M ; ++i ) {
      A( j , i ) = 1;
    }
  }

  //Deep copy host view to device views
  Kokkos::deep_copy(devy, y);
  Kokkos::deep_copy(devx, x);
  Kokkos::deep_copy(devA, A);

  // Timer products
  struct timeval begin,end;

  gettimeofday(&begin,NULL);

  for ( int repeat = 0; repeat < nrepeat; repeat++) {

    //Application: <y,Ax> = y^T*A*x
    double result = 0;
    Kokkos::parallel_reduce( range_policy( 0, N ), KOKKOS_LAMBDA ( int j, double &update ) {
      double temp2 = 0;
      for ( int i = 0 ; i < M ; ++i ) {
        temp2 += devA( j , i ) * devx( i );
      }
      update += devy( j ) * temp2;
    }, result );

    //Output result
    if ( repeat == (nrepeat - 1) )
      printf("  Computed result for %d x %d is %lf\n", N, M, result);
    const double solution = (double)N *(double)M;

    if ( result != solution ) {
      printf("  Error: result( %lf ) != solution( %lf )\n",result,solution);
    }
  }

  gettimeofday(&end,NULL);

  // Calculate time
  double time = 1.0*(end.tv_sec-begin.tv_sec) +
                1.0e-6*(end.tv_usec-begin.tv_usec);

  // Calculate bandwidth.
  // Each matrix A row (each of length M) is read once.
  // The x vector (of length M) is read N times.
  // The y vector (of length N) is read once.
  // double Gbytes = 1.0e-9 * double(sizeof(double) * ( 2 * M * N + N ));
  double Gbytes = 1.0e-9 * double(sizeof(double) * ( M + M * N + N )) ;

  // Print results (problem size, time and bandwidth in GB/s)
  printf("  M( %d ) N( %d ) nrepeat ( %d ) problem( %g MB ) time( %g s ) bandwidth( %g GB/s )\n",
         M , N, nrepeat, Gbytes * 1000, time, Gbytes * nrepeat / time );

  Kokkos::finalize();

  return 0 ;
}

However I get compiler errors with cuda. The error is the following.
[100%] Building CXX object CMakeFiles/trilinosExample.dir/ex1.cpp.o
In file included from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/impl/Kokkos_Atomic_View.hpp:47:0,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/impl/KokkosExp_ViewMapping.hpp:56,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_View.hpp:385,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Parallel.hpp:52,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Serial.hpp:52,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Core.hpp:53,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/cmake_cpp_examples/ex1.cpp:7:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Atomic.hpp:129:1: error: ‘__device__’ does not name a type
 __device__ inline
 ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Atomic.hpp:138:1: error: ‘__device__’ does not name a type
 __device__ inline
 ^
In file included from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Cuda.hpp:56:0,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Core.hpp:65,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/cmake_cpp_examples/ex1.cpp:7:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:574:10: error: ‘cudaTextureObject_t’ in namespace ‘::’ does not name a type
   static ::cudaTextureObject_t
          ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:581:3: error: ‘cudaTextureObject_t’ in namespace ‘::’ does not name a type
   ::cudaTextureObject_t   m_tex_obj ;
   ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:622:3: error: ‘cudaTextureObject_t’ in namespace ‘::’ does not name a type
   ::cudaTextureObject_t attach_texture_object()
   ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp: In constructor ‘Kokkos::Experimental::Impl::SharedAllocationRecord<Kokkos::CudaSpace, void>::SharedAllocationRecord()’:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:587:44: error: class ‘Kokkos::Experimental::Impl::SharedAllocationRecord<Kokkos::CudaSpace, void>’ does not have any field named ‘m_tex_obj’
   SharedAllocationRecord() : RecordBase(), m_tex_obj(0), m_space() {}
                                            ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp: At global scope:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:665:3: error: ‘cudaTextureObject_t’ in namespace ‘::’ does not name a type
   ::cudaTextureObject_t      m_tex_obj ;
   ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:708:3: error: ‘cudaTextureObject_t’ in namespace ‘::’ does not name a type
   ::cudaTextureObject_t attach_texture_object()
   ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp: In constructor ‘Kokkos::Experimental::Impl::SharedAllocationRecord<Kokkos::CudaUVMSpace, void>::SharedAllocationRecord()’:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_CudaSpace.hpp:671:44: error: class ‘Kokkos::Experimental::Impl::SharedAllocationRecord<Kokkos::CudaUVMSpace, void>’ does not have any field named ‘m_tex_obj’
   SharedAllocationRecord() : RecordBase(), m_tex_obj(0), m_space() {}
                                            ^
In file included from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Core.hpp:65:0,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/cmake_cpp_examples/ex1.cpp:7:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Cuda.hpp: At global scope:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Cuda.hpp:212:3: error: ‘cudaStream_t’ does not name a type
   cudaStream_t cuda_stream() const { return m_stream ; }
   ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Cuda.hpp:220:3: error: ‘cudaStream_t’ does not name a type
   cudaStream_t m_stream ;
   ^
In file included from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_CudaExec.hpp:56:0,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Cuda.hpp:260,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Kokkos_Core.hpp:65,
                 from /home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/cmake_cpp_examples/ex1.cpp:7:
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:56:33: error: variable or field ‘cuda_internal_error_throw’ declared void
 void cuda_internal_error_throw( cudaError e , const char * name, const char * file = NULL, const int line = 0 );
                                 ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:56:33: error: ‘cudaError’ was not declared in this scope
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:56:47: error: expected primary-expression before ‘const’
 void cuda_internal_error_throw( cudaError e , const char * name, const char * file = NULL, const int line = 0 );
                                               ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:56:66: error: expected primary-expression before ‘const’
 void cuda_internal_error_throw( cudaError e , const char * name, const char * file = NULL, const int line = 0 );
                                                                  ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:56:92: error: expected primary-expression before ‘const’
 void cuda_internal_error_throw( cudaError e , const char * name, const char * file = NULL, const int line = 0 );
                                                                                            ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:58:38: error: variable or field ‘cuda_internal_safe_call’ declared void
 inline void cuda_internal_safe_call( cudaError e , const char * name, const char * file = NULL, const int line = 0)
                                      ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:58:38: error: ‘cudaError’ was not declared in this scope
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:58:52: error: expected primary-expression before ‘const’
 inline void cuda_internal_safe_call( cudaError e , const char * name, const char * file = NULL, const int line = 0)
                                                    ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:58:71: error: expected primary-expression before ‘const’
 inline void cuda_internal_safe_call( cudaError e , const char * name, const char * file = NULL, const int line = 0)
                                                                       ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/install/include/Cuda/Kokkos_Cuda_Error.hpp:58:97: error: expected primary-expression before ‘const’
 inline void cuda_internal_safe_call( cudaError e , const char * name, const char * file = NULL, const int line = 0)
                                                                                                 ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/cmake_cpp_examples/ex1.cpp:124:1: error: expected ‘}’ at end of input
 }
 ^
/home/ashesh/Ashesh/Trilinos-Kokkos_GPU/Trilinos/cmake_cpp_examples/ex1.cpp:124:1: error: expected ‘}’ at end of input
make[2]: *** [CMakeFiles/trilinosExample.dir/ex1.cpp.o] Error 1
make[1]: *** [CMakeFiles/trilinosExample.dir/all] Error 2
make: *** [all] Error 2
It appears that there are some errors referenced from the hpp files. I want to know if I am making a mistake in linking with CUDA.
Tha CMakeList.txt file I use is attached.




--
Research Assistant, Multi-Scale Multi-Physics Computation Lab
University of Texas at El Paso
Texas, El Paso, USA
(+1) 915-355-5013
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20161101/388b378e/attachment.html>


More information about the Trilinos-Users mailing list