aboutsummaryrefslogtreecommitdiffhomepage
path: root/bench/btl/libs/gmm
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2008-08-04 23:12:48 +0000
committerGravatar Gael Guennebaud <g.gael@free.fr>2008-08-04 23:12:48 +0000
commita7a05382d1c51964bf3ea0536c6ddd9cc9888b72 (patch)
tree413a9ecf342bf59ff685495eff70f999ee7803ac /bench/btl/libs/gmm
parentc2f8ecf46683adcab0db05199ee2ebe15e6ada4a (diff)
Add a LU decomposition action in BTL and various cleaning in BTL. For instance
all per plot settings have been moved to a single file, go_mean now takes an optional second argument "tiny" to generate plots for tiny matrices, and output of comparison information wrt to previous benchs (if any).
Diffstat (limited to 'bench/btl/libs/gmm')
-rw-r--r--bench/btl/libs/gmm/gmm_interface.hh18
-rw-r--r--bench/btl/libs/gmm/main.cpp3
2 files changed, 15 insertions, 6 deletions
diff --git a/bench/btl/libs/gmm/gmm_interface.hh b/bench/btl/libs/gmm/gmm_interface.hh
index 93f827f9b..cfcf537ac 100644
--- a/bench/btl/libs/gmm/gmm_interface.hh
+++ b/bench/btl/libs/gmm/gmm_interface.hh
@@ -123,14 +123,20 @@ public :
gmm::lower_tri_solve(L, X, false);
}
- static inline void hessenberg(const gene_matrix & X, gene_matrix & C, int N){
- gmm::copy(X,C);
- gmm::Hessenberg_reduction(C,X,false);
+ static inline void lu_decomp(const gene_matrix & X, gene_matrix & R, int N){
+ gmm::copy(X,R);
+ std::vector<int> ipvt(N);
+ gmm::lu_factor(R, ipvt);
}
- static inline void tridiagonalization(const gene_matrix & X, gene_matrix & C, int N){
- gmm::copy(X,C);
- gmm::Householder_tridiagonalization(C,X,false);
+ static inline void hessenberg(const gene_matrix & X, gene_matrix & R, int N){
+ gmm::copy(X,R);
+ gmm::Hessenberg_reduction(R,X,false);
+ }
+
+ static inline void tridiagonalization(const gene_matrix & X, gene_matrix & R, int N){
+ gmm::copy(X,R);
+ gmm::Householder_tridiagonalization(R,X,false);
}
};
diff --git a/bench/btl/libs/gmm/main.cpp b/bench/btl/libs/gmm/main.cpp
index b7c36fc9a..d5aebcdbf 100644
--- a/bench/btl/libs/gmm/main.cpp
+++ b/bench/btl/libs/gmm/main.cpp
@@ -20,6 +20,7 @@
#include "bench.hh"
#include "basic_actions.hh"
#include "action_hessenberg.hh"
+#include "action_lu_decomp.hh"
BTL_MAIN;
@@ -39,6 +40,8 @@ int main()
bench<Action_trisolve<gmm_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
//bench<Action_lu_solve<blitz_LU_solve_interface<REAL_TYPE> > >(MIN_LU,MAX_LU,NB_POINT);
+ bench<Action_lu_decomp<gmm_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
+
bench<Action_hessenberg<gmm_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);
bench<Action_tridiagonalization<gmm_interface<REAL_TYPE> > >(MIN_MM,MAX_MM,NB_POINT);