aboutsummaryrefslogtreecommitdiffhomepage
path: root/bench/btl
diff options
context:
space:
mode:
Diffstat (limited to 'bench/btl')
-rw-r--r--bench/btl/CMakeLists.txt1
-rw-r--r--bench/btl/cmake/FindACML.cmake2
-rw-r--r--bench/btl/cmake/FindATLAS.cmake13
-rw-r--r--bench/btl/cmake/FindCBLAS.cmake1
-rw-r--r--bench/btl/cmake/FindOPENBLAS.cmake2
-rw-r--r--bench/btl/libs/eigen2/eigen2_interface.hh2
-rw-r--r--bench/btl/libs/tensors/CMakeLists.txt44
-rw-r--r--bench/btl/libs/tensors/main_linear.cpp23
-rw-r--r--bench/btl/libs/tensors/main_matmat.cpp21
-rw-r--r--bench/btl/libs/tensors/main_vecmat.cpp21
-rw-r--r--bench/btl/libs/tensors/tensor_interface.hh105
11 files changed, 224 insertions, 11 deletions
diff --git a/bench/btl/CMakeLists.txt b/bench/btl/CMakeLists.txt
index b299d9899..9444b450c 100644
--- a/bench/btl/CMakeLists.txt
+++ b/bench/btl/CMakeLists.txt
@@ -97,6 +97,7 @@ ENABLE_TESTING()
add_subdirectory(libs/eigen3)
add_subdirectory(libs/eigen2)
+add_subdirectory(libs/tensors)
add_subdirectory(libs/BLAS)
add_subdirectory(libs/ublas)
add_subdirectory(libs/gmm)
diff --git a/bench/btl/cmake/FindACML.cmake b/bench/btl/cmake/FindACML.cmake
index f45ae1b0d..4989fa2f4 100644
--- a/bench/btl/cmake/FindACML.cmake
+++ b/bench/btl/cmake/FindACML.cmake
@@ -17,6 +17,7 @@ find_file(ACML_LIBRARIES
libacml_mp.so
PATHS
/usr/lib
+ /usr/lib64
$ENV{ACMLDIR}/lib
${LIB_INSTALL_DIR}
)
@@ -35,6 +36,7 @@ if(NOT ACML_LIBRARIES)
libacml.so libacml_mv.so
PATHS
/usr/lib
+ /usr/lib64
$ENV{ACMLDIR}/lib
${LIB_INSTALL_DIR}
)
diff --git a/bench/btl/cmake/FindATLAS.cmake b/bench/btl/cmake/FindATLAS.cmake
index 14b1dee09..4136a989d 100644
--- a/bench/btl/cmake/FindATLAS.cmake
+++ b/bench/btl/cmake/FindATLAS.cmake
@@ -3,18 +3,13 @@ if (ATLAS_LIBRARIES)
set(ATLAS_FIND_QUIETLY TRUE)
endif (ATLAS_LIBRARIES)
-find_file(ATLAS_LIB libatlas.so.3 PATHS /usr/lib $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
+find_file(ATLAS_LIB libatlas.so.3 PATHS /usr/lib /usr/lib/atlas /usr/lib64 /usr/lib64/atlas $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
find_library(ATLAS_LIB satlas PATHS $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
-find_file(ATLAS_LAPACK liblapack_atlas.so.3 PATHS /usr/lib $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
-find_library(ATLAS_LAPACK lapack_atlas PATHS $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
+find_file(ATLAS_LAPACK NAMES liblapack_atlas.so.3 liblapack.so.3 PATHS /usr/lib /usr/lib/atlas /usr/lib64 /usr/lib64/atlas $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
+find_library(ATLAS_LAPACK NAMES lapack_atlas lapack PATHS $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
-if(NOT ATLAS_LAPACK)
- find_file(ATLAS_LAPACK liblapack.so.3 PATHS /usr/lib/atlas $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
- find_library(ATLAS_LAPACK lapack PATHS $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
-endif(NOT ATLAS_LAPACK)
-
-find_file(ATLAS_F77BLAS libf77blas.so.3 PATHS /usr/lib $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
+find_file(ATLAS_F77BLAS libf77blas.so.3 PATHS /usr/lib /usr/lib/atlas /usr/lib64 /usr/lib64/atlas $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
find_library(ATLAS_F77BLAS f77blas PATHS $ENV{ATLASDIR} ${LIB_INSTALL_DIR})
if(ATLAS_LIB AND ATLAS_CBLAS AND ATLAS_LAPACK AND ATLAS_F77BLAS)
diff --git a/bench/btl/cmake/FindCBLAS.cmake b/bench/btl/cmake/FindCBLAS.cmake
index 554f0291b..ce0f2f2b2 100644
--- a/bench/btl/cmake/FindCBLAS.cmake
+++ b/bench/btl/cmake/FindCBLAS.cmake
@@ -23,6 +23,7 @@ find_file(CBLAS_LIBRARIES
libcblas.so.3
PATHS
/usr/lib
+ /usr/lib64
$ENV{CBLASDIR}/lib
${LIB_INSTALL_DIR}
)
diff --git a/bench/btl/cmake/FindOPENBLAS.cmake b/bench/btl/cmake/FindOPENBLAS.cmake
index c76fc251c..2a0919436 100644
--- a/bench/btl/cmake/FindOPENBLAS.cmake
+++ b/bench/btl/cmake/FindOPENBLAS.cmake
@@ -3,7 +3,7 @@ if (OPENBLAS_LIBRARIES)
set(OPENBLAS_FIND_QUIETLY TRUE)
endif (OPENBLAS_LIBRARIES)
-find_file(OPENBLAS_LIBRARIES libopenblas.so PATHS /usr/lib $ENV{OPENBLASDIR} ${LIB_INSTALL_DIR})
+find_file(OPENBLAS_LIBRARIES NAMES libopenblas.so libopenblas.so.0 PATHS /usr/lib /usr/lib64 $ENV{OPENBLASDIR} ${LIB_INSTALL_DIR})
find_library(OPENBLAS_LIBRARIES openblas PATHS $ENV{OPENBLASDIR} ${LIB_INSTALL_DIR})
if(OPENBLAS_LIBRARIES AND CMAKE_COMPILER_IS_GNUCXX)
diff --git a/bench/btl/libs/eigen2/eigen2_interface.hh b/bench/btl/libs/eigen2/eigen2_interface.hh
index 47fe58135..1deabdae2 100644
--- a/bench/btl/libs/eigen2/eigen2_interface.hh
+++ b/bench/btl/libs/eigen2/eigen2_interface.hh
@@ -47,7 +47,7 @@ public :
{
#if defined(EIGEN_VECTORIZE_SSE)
if (SIZE==Dynamic) return "eigen2"; else return "tiny_eigen2";
- #elif defined(EIGEN_VECTORIZE_ALTIVEC)
+ #elif defined(EIGEN_VECTORIZE_ALTIVEC) || defined(EIGEN_VECTORIZE_VSX)
if (SIZE==Dynamic) return "eigen2"; else return "tiny_eigen2";
#else
if (SIZE==Dynamic) return "eigen2_novec"; else return "tiny_eigen2_novec";
diff --git a/bench/btl/libs/tensors/CMakeLists.txt b/bench/btl/libs/tensors/CMakeLists.txt
new file mode 100644
index 000000000..09d6d8e43
--- /dev/null
+++ b/bench/btl/libs/tensors/CMakeLists.txt
@@ -0,0 +1,44 @@
+
+
+if((NOT TENSOR_INCLUDE_DIR) AND Eigen_SOURCE_DIR)
+ # unless TENSOR_INCLUDE_DIR is defined, let's use current Eigen version
+ set(TENSOR_INCLUDE_DIR ${Eigen_SOURCE_DIR})
+ set(TENSOR_FOUND TRUE)
+else()
+ find_package(Tensor)
+endif()
+
+if (TENSOR_FOUND)
+
+ include_directories(${TENSOR_INCLUDE_DIR})
+ btl_add_bench(btl_tensor_linear main_linear.cpp)
+ btl_add_bench(btl_tensor_vecmat main_vecmat.cpp)
+ btl_add_bench(btl_tensor_matmat main_matmat.cpp)
+
+ btl_add_target_property(btl_tensor_linear COMPILE_FLAGS "-fno-exceptions -DBTL_PREFIX=tensor")
+ btl_add_target_property(btl_tensor_vecmat COMPILE_FLAGS "-fno-exceptions -DBTL_PREFIX=tensor")
+ btl_add_target_property(btl_tensor_matmat COMPILE_FLAGS "-fno-exceptions -DBTL_PREFIX=tensor")
+
+ option(BTL_BENCH_NOGCCVEC "also bench Eigen explicit vec without GCC's auto vec" OFF)
+ if(CMAKE_COMPILER_IS_GNUCXX AND BTL_BENCH_NOGCCVEC)
+ btl_add_bench(btl_tensor_nogccvec_linear main_linear.cpp)
+ btl_add_bench(btl_tensor_nogccvec_vecmat main_vecmat.cpp)
+ btl_add_bench(btl_tensor_nogccvec_matmat main_matmat.cpp)
+
+ btl_add_target_property(btl_tensor_nogccvec_linear COMPILE_FLAGS "-fno-exceptions -fno-tree-vectorize -DBTL_PREFIX=tensor_nogccvec")
+ btl_add_target_property(btl_tensor_nogccvec_vecmat COMPILE_FLAGS "-fno-exceptions -fno-tree-vectorize -DBTL_PREFIX=tensor_nogccvec")
+ btl_add_target_property(btl_tensor_nogccvec_matmat COMPILE_FLAGS "-fno-exceptions -fno-tree-vectorize -DBTL_PREFIX=tensor_nogccvec")
+ endif()
+
+
+ if(NOT BTL_NOVEC)
+ btl_add_bench(btl_tensor_novec_linear main_linear.cpp OFF)
+ btl_add_bench(btl_tensor_novec_vecmat main_vecmat.cpp OFF)
+ btl_add_bench(btl_tensor_novec_matmat main_matmat.cpp OFF)
+ btl_add_target_property(btl_tensor_novec_linear COMPILE_FLAGS "-fno-exceptions -DEIGEN_DONT_VECTORIZE -DBTL_PREFIX=tensor_novec")
+ btl_add_target_property(btl_tensor_novec_vecmat COMPILE_FLAGS "-fno-exceptions -DEIGEN_DONT_VECTORIZE -DBTL_PREFIX=tensor_novec")
+ btl_add_target_property(btl_tensor_novec_matmat COMPILE_FLAGS "-fno-exceptions -DEIGEN_DONT_VECTORIZE -DBTL_PREFIX=tensor_novec")
+
+ endif(NOT BTL_NOVEC)
+
+endif (TENSOR_FOUND)
diff --git a/bench/btl/libs/tensors/main_linear.cpp b/bench/btl/libs/tensors/main_linear.cpp
new file mode 100644
index 000000000..e257f1e72
--- /dev/null
+++ b/bench/btl/libs/tensors/main_linear.cpp
@@ -0,0 +1,23 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2014 Benoit Steiner <benoit.steiner.goog@gmail.com>
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+
+#include "utilities.h"
+#include "tensor_interface.hh"
+#include "bench.hh"
+#include "basic_actions.hh"
+
+BTL_MAIN;
+
+int main()
+{
+ bench<Action_axpy<tensor_interface<REAL_TYPE> > >(MIN_AXPY,MAX_AXPY,NB_POINT);
+ bench<Action_axpby<tensor_interface<REAL_TYPE> > >(MIN_AXPY,MAX_AXPY,NB_POINT);
+
+ return 0;
+}
diff --git a/bench/btl/libs/tensors/main_matmat.cpp b/bench/btl/libs/tensors/main_matmat.cpp
new file mode 100644
index 000000000..675fcfc6d
--- /dev/null
+++ b/bench/btl/libs/tensors/main_matmat.cpp
@@ -0,0 +1,21 @@
+//=====================================================
+// Copyright (C) 2014 Benoit Steiner <benoit.steiner.goog@gmail.com>
+//=====================================================
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+//
+#include "utilities.h"
+#include "tensor_interface.hh"
+#include "bench.hh"
+#include "basic_actions.hh"
+
+BTL_MAIN;
+
+int main()
+{
+ bench<Action_matrix_matrix_product<tensor_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
+
+ return 0;
+}
diff --git a/bench/btl/libs/tensors/main_vecmat.cpp b/bench/btl/libs/tensors/main_vecmat.cpp
new file mode 100644
index 000000000..1af00c81b
--- /dev/null
+++ b/bench/btl/libs/tensors/main_vecmat.cpp
@@ -0,0 +1,21 @@
+//=====================================================
+// Copyright (C) 2014 Benoit Steiner <benoit.steiner.goog@gmail.com>
+//=====================================================
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+//
+#include "utilities.h"
+#include "tensor_interface.hh"
+#include "bench.hh"
+#include "basic_actions.hh"
+
+BTL_MAIN;
+
+int main()
+{
+ bench<Action_matrix_vector_product<tensor_interface<REAL_TYPE> > >(MIN_MV,MAX_MV,NB_POINT);
+
+ return 0;
+}
diff --git a/bench/btl/libs/tensors/tensor_interface.hh b/bench/btl/libs/tensors/tensor_interface.hh
new file mode 100644
index 000000000..97b8e0f0b
--- /dev/null
+++ b/bench/btl/libs/tensors/tensor_interface.hh
@@ -0,0 +1,105 @@
+//=====================================================
+// Copyright (C) 2014 Benoit Steiner <benoit.steiner.goog@gmail.com>
+//=====================================================
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+//
+#ifndef TENSOR_INTERFACE_HH
+#define TENSOR_INTERFACE_HH
+
+#include <unsupported/Eigen/CXX11/Tensor>
+#include <vector>
+#include "btl.hh"
+
+using namespace Eigen;
+
+template<class real>
+class tensor_interface
+{
+public :
+ typedef real real_type;
+ typedef typename Eigen::Tensor<real,2>::Index Index;
+
+ typedef std::vector<real> stl_vector;
+ typedef std::vector<stl_vector> stl_matrix;
+
+ typedef Eigen::Tensor<real,2> gene_matrix;
+ typedef Eigen::Tensor<real,1> gene_vector;
+
+
+ static inline std::string name( void )
+ {
+ return EIGEN_MAKESTRING(BTL_PREFIX);
+ }
+
+ static void free_matrix(gene_matrix & /*A*/, int /*N*/) {}
+
+ static void free_vector(gene_vector & /*B*/) {}
+
+ static BTL_DONT_INLINE void matrix_from_stl(gene_matrix & A, stl_matrix & A_stl){
+ A.resize(Eigen::array<Index,2>(A_stl[0].size(), A_stl.size()));
+
+ for (unsigned int j=0; j<A_stl.size() ; j++){
+ for (unsigned int i=0; i<A_stl[j].size() ; i++){
+ A.coeffRef(Eigen::array<Index,2>(i,j)) = A_stl[j][i];
+ }
+ }
+ }
+
+ static BTL_DONT_INLINE void vector_from_stl(gene_vector & B, stl_vector & B_stl){
+ B.resize(B_stl.size());
+
+ for (unsigned int i=0; i<B_stl.size() ; i++){
+ B.coeffRef(i) = B_stl[i];
+ }
+ }
+
+ static BTL_DONT_INLINE void vector_to_stl(gene_vector & B, stl_vector & B_stl){
+ for (unsigned int i=0; i<B_stl.size() ; i++){
+ B_stl[i] = B.coeff(i);
+ }
+ }
+
+ static BTL_DONT_INLINE void matrix_to_stl(gene_matrix & A, stl_matrix & A_stl){
+ int N=A_stl.size();
+
+ for (int j=0;j<N;j++){
+ A_stl[j].resize(N);
+ for (int i=0;i<N;i++){
+ A_stl[j][i] = A.coeff(Eigen::array<Index,2>(i,j));
+ }
+ }
+ }
+
+ static inline void matrix_matrix_product(const gene_matrix & A, const gene_matrix & B, gene_matrix & X, int /*N*/){
+ typedef typename Eigen::Tensor<real_type, 1>::DimensionPair DimPair;
+ const Eigen::array<DimPair, 1> dims(DimPair(1, 0));
+ X/*.noalias()*/ = A.contract(B, dims);
+ }
+
+ static inline void matrix_vector_product(const gene_matrix & A, const gene_vector & B, gene_vector & X, int /*N*/){
+ typedef typename Eigen::Tensor<real_type, 1>::DimensionPair DimPair;
+ const Eigen::array<DimPair, 1> dims(DimPair(1, 0));
+ X/*.noalias()*/ = A.contract(B, dims);
+ }
+
+ static inline void axpy(real coef, const gene_vector & X, gene_vector & Y, int /*N*/){
+ Y += X.constant(coef) * X;
+ }
+
+ static inline void axpby(real a, const gene_vector & X, real b, gene_vector & Y, int /*N*/){
+ Y = X.constant(a)*X + Y.constant(b)*Y;
+ }
+
+ static EIGEN_DONT_INLINE void copy_matrix(const gene_matrix & source, gene_matrix & cible, int /*N*/){
+ cible = source;
+ }
+
+ static EIGEN_DONT_INLINE void copy_vector(const gene_vector & source, gene_vector & cible, int /*N*/){
+ cible = source;
+ }
+};
+
+#endif