Commit 1d894e61 authored by Tomáš Oberhuber's avatar Tomáš Oberhuber Committed by Jakub Klinkovský
Browse files

Fixed LightSpMV benchmark.

parent 7955b02a
Loading
Loading
Loading
Loading
+17 −3
Original line number Diff line number Diff line
@@ -61,17 +61,21 @@ struct LightSpMVBenchmark
   {
      static_assert( std::is_same< typename Matrix::DeviceType, TNL::Devices::Host >::value, "The only device type accepted here is TNL::Devices::Host." );
#ifdef HAVE_CUDA
      Options opt;
      cudaDeviceProp prop;
      cudaGetDeviceProperties(&prop, 0);
      opt._gpus.push_back(make_pair(0, prop));
      opt._numGPUs = 1;
      opt._numRows = matrix.getRows();
      opt._numCols = matrix.getColumns();
      opt._rowOffsets = matrix.getRowPointers().getData();
      opt._numValues = matrix.getValues().getSize();
      opt._colIndexValues = matrix.getColumnIndexes().getData();
      opt._numericalValues = matrix.getValues().getData();;
      opt._numericalValues = matrix.getValues().getData();
      opt._alpha = 1.0; // matrix multiplicator
      opt._beta = 0.0;  // output vector multiplicator
      opt._vectorX = inVector.getData();
      opt._vectorY = outVector.getData();
      opt._formula = 0;
      if( std::is_same< Real, float >::value )
      {
         if( kernelType == LightSpMVBenchmarkKernelVector )
@@ -130,7 +134,11 @@ struct LightSpMVBenchmark

   void vectorProduct()
   {
      this->spmv->invokeKernel( 0 );
#ifdef HAVE_CUDA
      this->spmv->spmvKernel();
      cudaDeviceSynchronize();
#endif

   }

   const CudaVectorView& getCudaOutVector()
@@ -142,11 +150,17 @@ struct LightSpMVBenchmark
   {
#ifdef HAVE_CUDA
      if( spmv ) delete spmv;
      opt._rowOffsets = nullptr;
      opt._colIndexValues = nullptr;
      opt._numericalValues = nullptr;
      opt._vectorX = nullptr;
      opt._vectorY = nullptr;
#endif
   }

   protected:
#ifdef HAVE_CUDA
      Options opt;
      SpMV* spmv = nullptr;
#endif
      VectorType  inVector, outVector;
+16 −5
Original line number Diff line number Diff line
@@ -38,9 +38,14 @@
#include <TNL/Algorithms/Segments/ChunkedEllpack.h>
#include <TNL/Algorithms/Segments/BiEllpack.h>

// Comment the following to turn off some groups of SpMV benchmarks and speed-up the compilation
#define WITH_TNL_BENCHMARK_SPMV_GENERAL_MATRICES
#define WITH_TNL_BENCHMARK_SPMV_SYMMETRIC_MATRICES
#define WITH_TNL_BENCHMARK_SPMV_LEGACY_FORMATS

// Uncomment the following line to enable benchmarking the sandbox sparse matrix.
//#define WITH_SANDBOX_MATRIX_BENCHMARK
#ifdef WITH_SANDBOX_MATRIX_BENCHMARK
//#define WITH_TNL_BENCHMARK_SPMV_SANDBOX_MATRIX
#ifdef WITH_TNL_BENCHMARK_SPMV_SANDBOX_MATRIX
#include <TNL/Matrices/Sandbox/SparseSandboxMatrix.h>
#endif

@@ -132,7 +137,7 @@ using BiEllpackSegments = Algorithms::Segments::BiEllpack< Device, Index, IndexA
template< typename Real, typename Device, typename Index >
using SymmetricSparseMatrix_BiEllpack = Matrices::SparseMatrix< Real, Device, Index, Matrices::SymmetricMatrix, BiEllpackSegments >;

#ifdef WITH_SANDBOX_MATRIX_BENCHMARK
#ifdef WITH_TNL_BENCHMARK_SPMV_SANDBOX_MATRIX
template< typename Real, typename Device, typename Index >
using SparseSandboxMatrix = Matrices::Sandbox::SparseSandboxMatrix< Real, Device, Index, Matrices::GeneralMatrix >;
#endif
@@ -491,6 +496,7 @@ benchmarkSpmvSynthetic( Benchmark& benchmark,
#endif
   csrHostMatrix.reset();

#ifdef WITH_TNL_BENCHMARK_SPMV_LEGACY_FORMATS
   /////
   // Benchmarking of TNL legacy formats
   //
@@ -515,7 +521,9 @@ benchmarkSpmvSynthetic( Benchmark& benchmark,
   }
   // AdEllpack is broken
   //benchmarkSpMV< Real, Matrices::AdEllpack              >( benchmark, hostOutVector, inputFileName, verboseMR );
#endif

#ifdef WITH_TNL_BENCHMARK_SPMV_GENERAL_MATRICES
   /////
   // Benchmarking TNL formats
   //
@@ -530,11 +538,13 @@ benchmarkSpmvSynthetic( Benchmark& benchmark,
   benchmarkSpMV< Real, HostMatrixType, SparseMatrix_SlicedEllpack                >( benchmark, hostMatrix, hostOutVector, inputFileName, verboseMR );
   benchmarkSpMV< Real, HostMatrixType, SparseMatrix_ChunkedEllpack               >( benchmark, hostMatrix, hostOutVector, inputFileName, verboseMR );
   benchmarkSpMV< Real, HostMatrixType, SparseMatrix_BiEllpack                    >( benchmark, hostMatrix, hostOutVector, inputFileName, verboseMR );
#ifdef WITH_SANDBOX_MATRIX_BENCHMARK
#ifdef WITH_TNL_BENCHMARK_SPMV_SANDBOX_MATRIX
   benchmarkSpMV< Real, HostMatrixType, SparseSandboxMatrix                       >( benchmark, hostMatrix, hostOutVector, inputFileName, verboseMR );
#endif
   hostMatrix.reset();
#endif

#ifdef WITH_TNL_BENCHMARK_SPMV_SYMMETRIC_MATRICES
   /////
   // Benchmarking symmetric sparse matrices
   //
@@ -556,7 +566,7 @@ benchmarkSpmvSynthetic( Benchmark& benchmark,
      TNL::Matrices::MatrixReader< InputMatrix >::readMtx( inputFileName, hostMatrix, verboseMR );
      if( hostMatrix != symmetricHostMatrix )
      {
         std::cerr << "ERROR !!!!!! " << std::endl;
         std::cerr << "ERROR: Symmetric matrices do not match !!!" << std::endl;
      }
      benchmarkSpMV< Real, SymmetricInputMatrix, SymmetricSparseMatrix_CSR_Scalar                   >( benchmark, symmetricHostMatrix, hostOutVector, inputFileName, verboseMR );
      benchmarkSpMV< Real, SymmetricInputMatrix, SymmetricSparseMatrix_CSR_Vector                   >( benchmark, symmetricHostMatrix, hostOutVector, inputFileName, verboseMR );
@@ -567,6 +577,7 @@ benchmarkSpmvSynthetic( Benchmark& benchmark,
      benchmarkSpMV< Real, SymmetricInputMatrix, SymmetricSparseMatrix_ChunkedEllpack               >( benchmark, symmetricHostMatrix, hostOutVector, inputFileName, verboseMR );
      benchmarkSpMV< Real, SymmetricInputMatrix, SymmetricSparseMatrix_BiEllpack                    >( benchmark, symmetricHostMatrix, hostOutVector, inputFileName, verboseMR );
   }
#endif
}

} // namespace SpMVLegacy
+9 −9
Original line number Diff line number Diff line
@@ -66,7 +66,7 @@ do
	 mtx_file_name=${mtx_file_name%.mtx}
         if test x$DEBUG = xyes;
         then
            gdb --args $BENCHMARK --input-file $file --log-file log-files/sparse-matrix-benchmark.log --output-mode append --verbose 1
            gdb --args ${BENCHMARK_DBG} --input-file $file --log-file log-files/sparse-matrix-benchmark.log --output-mode append --verbose 1
         else
            $BENCHMARK --input-file $file --log-file log-files/sparse-matrix-benchmark.log --output-mode append --verbose 1
         fi