summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorAndrea Arteaga <andyspiros@gmail.com>2012-10-07 23:36:29 +0200
committerAndrea Arteaga <andyspiros@gmail.com>2012-10-07 23:36:29 +0200
commit9dddb5af4eb913b6b607a773508f9ebb62c28c46 (patch)
treecf53caa402129ff1df5e16838a9a5bd065d088a7
parentAdded and tested actions for BLAS level 3. (diff)
downloadauto-numerical-bench-9dddb5af4eb913b6b607a773508f9ebb62c28c46.tar.gz
auto-numerical-bench-9dddb5af4eb913b6b607a773508f9ebb62c28c46.tar.bz2
auto-numerical-bench-9dddb5af4eb913b6b607a773508f9ebb62c28c46.zip
Placed BLAS action under subfolder. Written first LAPACK action. Updated
samples for blas and cblas.
-rw-r--r--btl/NumericInterface/NI_internal/FortranDeclarations.hpp14
-rw-r--r--btl/NumericInterface/NI_internal/FortranInterface.hpp12
-rw-r--r--btl/actions/BLAS/action_MatrixMatrix.hpp (renamed from btl/actions/action_MatrixMatrix.hpp)0
-rw-r--r--btl/actions/BLAS/action_MatrixTMatrix.hpp (renamed from btl/actions/action_MatrixTMatrix.hpp)0
-rw-r--r--btl/actions/BLAS/action_MatrixTVector.hpp (renamed from btl/actions/action_MatrixTVector.hpp)0
-rw-r--r--btl/actions/BLAS/action_MatrixVector.hpp (renamed from btl/actions/action_MatrixVector.hpp)0
-rw-r--r--btl/actions/BLAS/action_Rank1Update.hpp (renamed from btl/actions/action_Rank1Update.hpp)0
-rw-r--r--btl/actions/BLAS/action_Rank2Update.hpp (renamed from btl/actions/action_Rank2Update.hpp)0
-rw-r--r--btl/actions/BLAS/action_SymMatrixVector.hpp (renamed from btl/actions/action_SymMatrixVector.hpp)0
-rw-r--r--btl/actions/BLAS/action_TriMatrixMatrix.hpp (renamed from btl/actions/action_TriMatrixMatrix.hpp)0
-rw-r--r--btl/actions/BLAS/action_TriSolveMatrix.hpp (renamed from btl/actions/action_TriSolveMatrix.hpp)0
-rw-r--r--btl/actions/BLAS/action_TriSolveVector.hpp (renamed from btl/actions/action_TriSolveVector.hpp)0
-rw-r--r--btl/actions/BLAS/action_axpy.hpp (renamed from btl/actions/action_axpy.hpp)0
-rw-r--r--btl/actions/BLAS/action_rot.hpp (renamed from btl/actions/action_rot.hpp)0
-rw-r--r--btl/actions/LAPACK/GeneralSolve.hpp84
-rw-r--r--btl/actions/actionsBLAS.hpp41
-rw-r--r--samples/blastests.xml33
-rw-r--r--samples/cblastests.xml12
18 files changed, 151 insertions, 45 deletions
diff --git a/btl/NumericInterface/NI_internal/FortranDeclarations.hpp b/btl/NumericInterface/NI_internal/FortranDeclarations.hpp
index 35bce39..2378cf2 100644
--- a/btl/NumericInterface/NI_internal/FortranDeclarations.hpp
+++ b/btl/NumericInterface/NI_internal/FortranDeclarations.hpp
@@ -79,6 +79,20 @@ extern "C" {
void strsm_(const char*, const char*, const char*, const char*, const int*, const int*, const float*, const float*, const int*, float*, const int*);
void dtrsm_(const char*, const char*, const char*, const char*, const int*, const int*, const double*, const double*, const int*, double*, const int*);
+
+
+
+ /******************
+ * LAPACK SOLVERS *
+ ******************/
+
+ void sgesv(const int&, const int&, float*, const int&, int*, float*, const int&, int&);
+ void dgesv(const int&, const int&, double*, const int&, int*, double*, const int&, int&);
+
+
+
+ FORTFUNC(gesv)(&N, &ONE, A, &N, ipiv, b, &N, &info);
+
}
diff --git a/btl/NumericInterface/NI_internal/FortranInterface.hpp b/btl/NumericInterface/NI_internal/FortranInterface.hpp
index eea5898..2073e80 100644
--- a/btl/NumericInterface/NI_internal/FortranInterface.hpp
+++ b/btl/NumericInterface/NI_internal/FortranInterface.hpp
@@ -161,6 +161,18 @@ public:
{
FORTFUNC(trsm)("L", &uplo, "N", "N", &M, &N, &fONE, A, &M, B, &M);
}
+
+
+
+ /******************
+ * LAPACK SOLVERS *
+ ******************/
+
+ static void GeneralSolve(const int& N, Scalar *A, Scalar *b, int *ipiv)
+ {
+ int info;
+ FORTFUNC(gesv)(&N, &ONE, A, &N, ipiv, b, &N, &info);
+ }
};
const int NumericInterface<NI_SCALAR>::ZERO = 0;
diff --git a/btl/actions/action_MatrixMatrix.hpp b/btl/actions/BLAS/action_MatrixMatrix.hpp
index a69a7f0..a69a7f0 100644
--- a/btl/actions/action_MatrixMatrix.hpp
+++ b/btl/actions/BLAS/action_MatrixMatrix.hpp
diff --git a/btl/actions/action_MatrixTMatrix.hpp b/btl/actions/BLAS/action_MatrixTMatrix.hpp
index 43aaf74..43aaf74 100644
--- a/btl/actions/action_MatrixTMatrix.hpp
+++ b/btl/actions/BLAS/action_MatrixTMatrix.hpp
diff --git a/btl/actions/action_MatrixTVector.hpp b/btl/actions/BLAS/action_MatrixTVector.hpp
index 20b360a..20b360a 100644
--- a/btl/actions/action_MatrixTVector.hpp
+++ b/btl/actions/BLAS/action_MatrixTVector.hpp
diff --git a/btl/actions/action_MatrixVector.hpp b/btl/actions/BLAS/action_MatrixVector.hpp
index b1a37d0..b1a37d0 100644
--- a/btl/actions/action_MatrixVector.hpp
+++ b/btl/actions/BLAS/action_MatrixVector.hpp
diff --git a/btl/actions/action_Rank1Update.hpp b/btl/actions/BLAS/action_Rank1Update.hpp
index bd88ac3..bd88ac3 100644
--- a/btl/actions/action_Rank1Update.hpp
+++ b/btl/actions/BLAS/action_Rank1Update.hpp
diff --git a/btl/actions/action_Rank2Update.hpp b/btl/actions/BLAS/action_Rank2Update.hpp
index b0b1693..b0b1693 100644
--- a/btl/actions/action_Rank2Update.hpp
+++ b/btl/actions/BLAS/action_Rank2Update.hpp
diff --git a/btl/actions/action_SymMatrixVector.hpp b/btl/actions/BLAS/action_SymMatrixVector.hpp
index e475007..e475007 100644
--- a/btl/actions/action_SymMatrixVector.hpp
+++ b/btl/actions/BLAS/action_SymMatrixVector.hpp
diff --git a/btl/actions/action_TriMatrixMatrix.hpp b/btl/actions/BLAS/action_TriMatrixMatrix.hpp
index 2b8a7c1..2b8a7c1 100644
--- a/btl/actions/action_TriMatrixMatrix.hpp
+++ b/btl/actions/BLAS/action_TriMatrixMatrix.hpp
diff --git a/btl/actions/action_TriSolveMatrix.hpp b/btl/actions/BLAS/action_TriSolveMatrix.hpp
index f2f8879..f2f8879 100644
--- a/btl/actions/action_TriSolveMatrix.hpp
+++ b/btl/actions/BLAS/action_TriSolveMatrix.hpp
diff --git a/btl/actions/action_TriSolveVector.hpp b/btl/actions/BLAS/action_TriSolveVector.hpp
index 6cac6f1..6cac6f1 100644
--- a/btl/actions/action_TriSolveVector.hpp
+++ b/btl/actions/BLAS/action_TriSolveVector.hpp
diff --git a/btl/actions/action_axpy.hpp b/btl/actions/BLAS/action_axpy.hpp
index 85cb40e..85cb40e 100644
--- a/btl/actions/action_axpy.hpp
+++ b/btl/actions/BLAS/action_axpy.hpp
diff --git a/btl/actions/action_rot.hpp b/btl/actions/BLAS/action_rot.hpp
index 574254a..574254a 100644
--- a/btl/actions/action_rot.hpp
+++ b/btl/actions/BLAS/action_rot.hpp
diff --git a/btl/actions/LAPACK/GeneralSolve.hpp b/btl/actions/LAPACK/GeneralSolve.hpp
new file mode 100644
index 0000000..7fb98ba
--- /dev/null
+++ b/btl/actions/LAPACK/GeneralSolve.hpp
@@ -0,0 +1,84 @@
+//=====================================================
+// Copyright (C) 2012 Andrea Arteaga <andyspiros@gmail.com>
+//=====================================================
+//
+// This program is free software; you can redistribute it and/or
+// modify it under the terms of the GNU General Public License
+// as published by the Free Software Foundation; either version 2
+// of the License, or (at your option) any later version.
+//
+// This program is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+// GNU General Public License for more details.
+// You should have received a copy of the GNU General Public License
+// along with this program; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
+//
+#ifndef ACTION_GENERALSOLVE
+#define ACTION_GENERALSOLVE
+
+#include "LinearCongruential.hpp"
+#include <vector>
+#include <algorithm>
+
+template<class Interface>
+class Action_GeneralSolve {
+
+ typedef typename Interface::Scalar Scalar;
+ typedef std::vector<Scalar> vector_t;
+
+private:
+ // Invalidate copy constructor
+ Action_GeneralSolve(const Action_GeneralSolve&);
+
+public:
+
+ // Constructor
+ Action_GeneralSolve(int size)
+ : _size(size), lc(10),
+ A(lc.fillVector<Scalar>(size*size)), b(lc.fillVector<Scalar>(size)),
+ A_work(size*size), x_work(size), b_res(size), ipiv(size)
+ {
+ MESSAGE("Action_GeneralSolve Constructor");
+ }
+
+ // Action name
+ static std::string name()
+ {
+ return "GeneralSolve_" + Interface::name();
+ }
+
+ double fpo() {
+ return double(_size)*double(_size)*double(_size);
+ }
+
+ inline void initialize(){
+ std::copy(A.begin(), A.end(), A_work.begin());
+ std::copy(b.begin(), b.end(), x_work.begin());
+ }
+
+ inline void calculate() {
+ Interface::GeneralSolve(_size, &A_work[0], &x_work[0], &ipiv[0]);
+ }
+
+ Scalar getResidual() {
+ initialize();
+ calculate();
+ std::copy(b.begin(), b.end(), b_res.begin());
+ Interface::MatrixVector(false, _size, _size, -1., &A[0], &x_work[0],
+ 1., &b_res[0]);
+ return Interface::norm(_size, &b_res[0]);
+ }
+
+private:
+ const int _size;
+ LinearCongruential<> lc;
+
+ const vector_t A, b;
+ vector_t A_work, x_work, b_res;
+ std::vector<int> ipiv;
+
+};
+
+#endif // ACTION_GENERALSOLVE
diff --git a/btl/actions/actionsBLAS.hpp b/btl/actions/actionsBLAS.hpp
index 1ac4b12..9d8aed8 100644
--- a/btl/actions/actionsBLAS.hpp
+++ b/btl/actions/actionsBLAS.hpp
@@ -1,14 +1,31 @@
-#include "action_axpy.hpp"
-#include "action_rot.hpp"
+//=====================================================
+// Copyright (C) 2012 Andrea Arteaga <andyspiros@gmail.com>
+//=====================================================
+//
+// This program is free software; you can redistribute it and/or
+// modify it under the terms of the GNU General Public License
+// as published by the Free Software Foundation; either version 2
+// of the License, or (at your option) any later version.
+//
+// This program is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+// GNU General Public License for more details.
+// You should have received a copy of the GNU General Public License
+// along with this program; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
+//
+#include "BLAS/action_axpy.hpp"
+#include "BLAS/action_rot.hpp"
-#include "action_MatrixVector.hpp"
-#include "action_MatrixTVector.hpp"
-#include "action_SymMatrixVector.hpp"
-#include "action_TriSolveVector.hpp"
-#include "action_Rank1Update.hpp"
-#include "action_Rank2Update.hpp"
+#include "BLAS/action_MatrixVector.hpp"
+#include "BLAS/action_MatrixTVector.hpp"
+#include "BLAS/action_SymMatrixVector.hpp"
+#include "BLAS/action_TriSolveVector.hpp"
+#include "BLAS/action_Rank1Update.hpp"
+#include "BLAS/action_Rank2Update.hpp"
-#include "action_MatrixMatrix.hpp"
-#include "action_MatrixTMatrix.hpp"
-#include "action_TriMatrixMatrix.hpp"
-#include "action_TriSolveMatrix.hpp"
+#include "BLAS/action_MatrixMatrix.hpp"
+#include "BLAS/action_MatrixTMatrix.hpp"
+#include "BLAS/action_TriMatrixMatrix.hpp"
+#include "BLAS/action_TriSolveMatrix.hpp"
diff --git a/samples/blastests.xml b/samples/blastests.xml
index c9f45c2..fba38f2 100644
--- a/samples/blastests.xml
+++ b/samples/blastests.xml
@@ -1,44 +1,23 @@
<numbench>
- <operations module="blas">axpy matrix_vector matrix_matrix aat</operations>
+ <operations module="blas">axpy MatrixVector MatrixTVector MatrixMatrix MatrixTMatrix TriSolveVector TriSolveMatrix</operations>
<testcases>
- <case id="reference">
+ <case id="reference-Os">
<pkg>sci-libs/blas-reference-20120423</pkg>
<emergeenv>
- <var name="FFLAGS">-O3</var>
+ <var name="FFLAGS">-Os</var>
</emergeenv>
</case>
- <case id="atlas">
- <pkg>sci-libs/atlas-3.10.0</pkg>
- <emergeenv>
- <var name="USE">threads</var>
- </emergeenv>
- </case>
-
- <case id="eigen">
- <pkg>dev-cpp/eigen-3.1.1-r1</pkg>
- <emergeenv>
- <var name="CXXFLAGS">-O3</var>
- </emergeenv>
- </case>
-
- <case id="openblas">
- <pkg>sci-libs/openblas-0.2.3</pkg>
+ <case id="reference-Ofast">
+ <pkg>sci-libs/blas-reference-20120423</pkg>
<emergeenv>
- <var name="FFLAGS">-O3</var>
- <var name="CFLAGS">-O3</var>
- <var name="USE">-openmp threads</var>
+ <var name="FFLAGS">-Ofast</var>
</emergeenv>
</case>
- <case id="mkl">
- <pkg>sci-libs/mkl-10.3.7.256</pkg>
- <skip>mkl32*</skip>
- </case>
-
</testcases>
</numbench>
diff --git a/samples/cblastests.xml b/samples/cblastests.xml
index 9a08eef..00c6644 100644
--- a/samples/cblastests.xml
+++ b/samples/cblastests.xml
@@ -1,13 +1,13 @@
<numbench>
- <operations module="cblas">axpy matrix_vector matrix_matrix aat</operations>
+ <operations module="cblas">axpy MatrixVector MatrixTVector MatrixMatrix MatrixTMatrix TriSolveVector TriSolveMatrix</operations>
<testcases>
<case id="reference">
<pkg>sci-libs/cblas-reference-20110218</pkg>
<emergeenv>
- <var name="FFLAGS">-O3</var>
+ <var name="FFLAGS">-Ofast</var>
</emergeenv>
</case>
@@ -18,22 +18,22 @@
</emergeenv>
</case>
- <case id="openblas1">
+ <case id="openblas">
<pkg>sci-libs/openblas-0.2.3</pkg>
<emergeenv>
<var name="FFLAGS">-O3</var>
<var name="CFLAGS">-O3</var>
<var name="USE">incblas</var>
</emergeenv>
- <runenv>
- <var name="OPENBLAS_NUM_THREADS">1</var>
- </runenv>
</case>
<case id="mkl">
<pkg>sci-libs/mkl-10.3.7.256</pkg>
<skip>mkl32*</skip>
<skip>*int64*</skip>
+ <skip>mkl64-dynamic</skip>
+ <skip>mkl64-gfortran</skip>
+ <skip>mkl64-intel</skip>
</case>
</testcases>