aboutsummaryrefslogtreecommitdiffhomepage
path: root/bench/btl/libs
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-06-04 18:16:54 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-06-04 18:16:54 +0200
commitc49d1fd2b5d3ef6a71c385a84755df94434eb4b3 (patch)
tree65726702657bc67caff17ecc814cb805b0f25bc9 /bench/btl/libs
parentf26c691678a7c6dd02789af3563a042e82677849 (diff)
add a partial LU bench in BTL
Diffstat (limited to 'bench/btl/libs')
-rw-r--r--bench/btl/libs/C_BLAS/C_BLAS_interface.hh12
-rw-r--r--bench/btl/libs/C_BLAS/main.cpp2
-rw-r--r--bench/btl/libs/eigen2/eigen2_interface.hh5
-rw-r--r--bench/btl/libs/eigen2/main_adv.cpp2
4 files changed, 19 insertions, 2 deletions
diff --git a/bench/btl/libs/C_BLAS/C_BLAS_interface.hh b/bench/btl/libs/C_BLAS/C_BLAS_interface.hh
index db02c2e2e..c85c4e818 100644
--- a/bench/btl/libs/C_BLAS/C_BLAS_interface.hh
+++ b/bench/btl/libs/C_BLAS/C_BLAS_interface.hh
@@ -61,6 +61,7 @@ extern "C"
void sgehrd_( const int *n, int *ilo, int *ihi, float *a, const int *lda, float *tau, float *work, int *lwork, int *info );
// LU row pivoting
+// void dgetrf_( int *m, int *n, double *a, int *lda, int *ipiv, int *info );
// void sgetrf_(const int* m, const int* n, float *a, const int* ld, int* ipivot, int* info);
// LU full pivoting
void sgetc2_(const int* n, float *a, const int *lda, int *ipiv, int *jpiv, int*info );
@@ -160,7 +161,7 @@ public :
cblas_ssymv(CblasColMajor,CblasLower,N,1.0,A,N,B,1,0.0,X,1);
#endif
}
-
+
static inline void syr2(gene_matrix & A, gene_vector & B, gene_vector & X, int N){
#ifdef PUREBLAS
ssyr2_(&lower,&N,&fone,B,&intone,X,&intone,A,&N);
@@ -247,6 +248,15 @@ public :
sgetc2_(&N, C, &N, ipiv, jpiv, &info);
}
+ static inline void partial_lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
+ int N2 = N*N;
+ scopy_(&N2, X, &intone, C, &intone);
+ char uplo = 'L';
+ int info = 0;
+ int * ipiv = (int*)alloca(sizeof(int)*N);
+ sgetrf_(&N, &N, C, &N, ipiv, &info);
+ }
+
static inline void hessenberg(const gene_matrix & X, gene_matrix & C, int N){
#ifdef PUREBLAS
{
diff --git a/bench/btl/libs/C_BLAS/main.cpp b/bench/btl/libs/C_BLAS/main.cpp
index 57cb9930e..81cac6c29 100644
--- a/bench/btl/libs/C_BLAS/main.cpp
+++ b/bench/btl/libs/C_BLAS/main.cpp
@@ -24,6 +24,7 @@
#include "action_cholesky.hh"
#include "action_lu_decomp.hh"
+#include "action_partial_lu_decomp.hh"
#include "action_trisolve_matrix.hh"
#ifdef HAS_LAPACK
@@ -54,6 +55,7 @@ int main()
#ifdef HAS_LAPACK
bench<Action_lu_decomp<C_BLAS_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
+ bench<Action_partial_lu_decomp<C_BLAS_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_hessenberg<C_BLAS_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_tridiagonalization<C_BLAS_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
#endif
diff --git a/bench/btl/libs/eigen2/eigen2_interface.hh b/bench/btl/libs/eigen2/eigen2_interface.hh
index 5d70279ca..186a485cf 100644
--- a/bench/btl/libs/eigen2/eigen2_interface.hh
+++ b/bench/btl/libs/eigen2/eigen2_interface.hh
@@ -202,7 +202,10 @@ public :
static inline void lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
C = X.lu().matrixLU();
-// C = X.inverse();
+ }
+
+ static inline void partial_lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
+ C = X.partialLu().matrixLU();
}
static inline void tridiagonalization(const gene_matrix & X, gene_matrix & C, int N){
diff --git a/bench/btl/libs/eigen2/main_adv.cpp b/bench/btl/libs/eigen2/main_adv.cpp
index 6dd9d3884..e6c7fe7cd 100644
--- a/bench/btl/libs/eigen2/main_adv.cpp
+++ b/bench/btl/libs/eigen2/main_adv.cpp
@@ -23,6 +23,7 @@
#include "action_cholesky.hh"
#include "action_hessenberg.hh"
#include "action_lu_decomp.hh"
+#include "action_partial_lu_decomp.hh"
BTL_MAIN;
@@ -32,6 +33,7 @@ int main()
bench<Action_trisolve_matrix<eigen2_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_cholesky<eigen2_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_lu_decomp<eigen2_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
+ bench<Action_partial_lu_decomp<eigen2_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_hessenberg<eigen2_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_tridiagonalization<eigen2_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);