summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorAlberto Ruiz <aruiz@um.es>2007-10-02 18:59:50 +0000
committerAlberto Ruiz <aruiz@um.es>2007-10-02 18:59:50 +0000
commit42bec1ac9911131b552f66779203eb599a86563d (patch)
treec4aefaedb21730644fcd4d2f85d830fe4d4daf07
parentd925bada507562250a75587c409bdb35bbbc6ed8 (diff)
lapack real and complex unpacked QR
-rw-r--r--examples/tests.hs14
-rw-r--r--examples/useStatic.hs (renamed from examples/usaStatic.hs)0
-rw-r--r--lib/Numeric/GSL/Matrix.hs26
-rw-r--r--lib/Numeric/GSL/Vector.hs1
-rw-r--r--lib/Numeric/GSL/gsl-aux.c27
-rw-r--r--lib/Numeric/GSL/gsl-aux.h3
-rw-r--r--lib/Numeric/LinearAlgebra/Algorithms.hs38
-rw-r--r--lib/Numeric/LinearAlgebra/LAPACK/lapack-aux.c5
8 files changed, 101 insertions, 13 deletions
diff --git a/examples/tests.hs b/examples/tests.hs
index b088069..b3272cf 100644
--- a/examples/tests.hs
+++ b/examples/tests.hs
@@ -7,7 +7,8 @@ import Numeric.GSL hiding (sin,cos,exp,choose)
7import Numeric.LinearAlgebra 7import Numeric.LinearAlgebra
8import Numeric.LinearAlgebra.Linear(Linear) 8import Numeric.LinearAlgebra.Linear(Linear)
9import Numeric.LinearAlgebra.LAPACK 9import Numeric.LinearAlgebra.LAPACK
10import Numeric.GSL.Matrix 10import Numeric.GSL.Matrix(svdg)
11import qualified Numeric.GSL.Matrix as GSL
11import Test.QuickCheck hiding (test) 12import Test.QuickCheck hiding (test)
12import Test.HUnit hiding ((~:),test) 13import Test.HUnit hiding ((~:),test)
13import System.Random(randomRs,mkStdGen) 14import System.Random(randomRs,mkStdGen)
@@ -271,6 +272,11 @@ cholCTest = chol ((2><2) [1,2,2,9::Complex Double]) == (2><2) [1,2,0,2.236067977
271 272
272--------------------------------------------------------------------- 273---------------------------------------------------------------------
273 274
275qrTest qr m = q <> r |~| m -- && q <> ctrans q |~| ident (rows m)
276 where (q,r) = qr m
277
278---------------------------------------------------------------------
279
274asFortran m = (rows m >|< cols m) $ toList (fdat m) 280asFortran m = (rows m >|< cols m) $ toList (fdat m)
275asC m = (rows m >< cols m) $ toList (cdat m) 281asC m = (rows m >< cols m) $ toList (cdat m)
276 282
@@ -322,6 +328,12 @@ tests = do
322 [ test "cholR" cholRTest 328 [ test "cholR" cholRTest
323 , test "cholC" cholRTest 329 , test "cholC" cholRTest
324 ] 330 ]
331 putStrLn "--------- qr ---------"
332 quickCheck (qrTest GSL.qr)
333 quickCheck (qrTest (GSL.unpackQR . GSL.qrPacked))
334 quickCheck (qrTest ( unpackQR . GSL.qrPacked))
335 quickCheck (qrTest qr ::RM->Bool)
336 quickCheck (qrTest qr ::CM->Bool)
325 putStrLn "--------- nullspace ------" 337 putStrLn "--------- nullspace ------"
326 quickCheck (nullspaceTest :: RM -> Bool) 338 quickCheck (nullspaceTest :: RM -> Bool)
327 quickCheck (nullspaceTest :: CM -> Bool) 339 quickCheck (nullspaceTest :: CM -> Bool)
diff --git a/examples/usaStatic.hs b/examples/useStatic.hs
index 619af8f..619af8f 100644
--- a/examples/usaStatic.hs
+++ b/examples/useStatic.hs
diff --git a/lib/Numeric/GSL/Matrix.hs b/lib/Numeric/GSL/Matrix.hs
index eb1931a..5a5c19e 100644
--- a/lib/Numeric/GSL/Matrix.hs
+++ b/lib/Numeric/GSL/Matrix.hs
@@ -16,7 +16,7 @@
16module Numeric.GSL.Matrix( 16module Numeric.GSL.Matrix(
17 eigSg, eigHg, 17 eigSg, eigHg,
18 svdg, 18 svdg,
19 qr, 19 qr, qrPacked, unpackQR,
20 cholR, -- cholC, 20 cholR, -- cholC,
21 luSolveR, luSolveC, 21 luSolveR, luSolveC,
22 luR, luC 22 luR, luC
@@ -149,6 +149,30 @@ qr x = unsafePerformIO $ do
149 c = cols x 149 c = cols x
150foreign import ccall "gsl-aux.h QR" c_qr :: TMMM 150foreign import ccall "gsl-aux.h QR" c_qr :: TMMM
151 151
152qrPacked :: Matrix Double -> (Matrix Double, Vector Double)
153qrPacked x = unsafePerformIO $ do
154 qr <- createMatrix RowMajor r c
155 tau <- createVector (min r c)
156 c_qrPacked // mat cdat x // mat dat qr // vec tau // check "qrUnpacked" [cdat x]
157 return (qr,tau)
158 where r = rows x
159 c = cols x
160foreign import ccall "gsl-aux.h QRpacked" c_qrPacked :: TMMV
161
162unpackQR :: (Matrix Double, Vector Double) -> (Matrix Double, Matrix Double)
163unpackQR (qr,tau) = unsafePerformIO $ do
164 q <- createMatrix RowMajor r r
165 rot <- createMatrix RowMajor r c
166 c_qrUnpack // mat cdat qr // vec tau // mat dat q // mat dat rot // check "qrUnpack" [cdat qr,tau]
167 return (q,rot)
168 where r = rows qr
169 c = cols qr
170foreign import ccall "gsl-aux.h QRunpack" c_qrUnpack :: TMVMM
171
172
173type TMMV = Int -> Int -> PD -> TMV
174type TMVMM = Int -> Int -> PD -> Int -> PD -> TMM
175
152{- | Cholesky decomposition of a symmetric positive definite real matrix using /gsl_linalg_cholesky_decomp/. 176{- | Cholesky decomposition of a symmetric positive definite real matrix using /gsl_linalg_cholesky_decomp/.
153 177
154@\> chol $ (2><2) [1,2, 178@\> chol $ (2><2) [1,2,
diff --git a/lib/Numeric/GSL/Vector.hs b/lib/Numeric/GSL/Vector.hs
index ef3d5e8..41efdc0 100644
--- a/lib/Numeric/GSL/Vector.hs
+++ b/lib/Numeric/GSL/Vector.hs
@@ -12,7 +12,6 @@
12-- Vector operations 12-- Vector operations
13-- 13--
14----------------------------------------------------------------------------- 14-----------------------------------------------------------------------------
15-- #hide
16 15
17module Numeric.GSL.Vector ( 16module Numeric.GSL.Vector (
18 FunCodeS(..), toScalarR, 17 FunCodeS(..), toScalarR,
diff --git a/lib/Numeric/GSL/gsl-aux.c b/lib/Numeric/GSL/gsl-aux.c
index c602d5e..83c3bf8 100644
--- a/lib/Numeric/GSL/gsl-aux.c
+++ b/lib/Numeric/GSL/gsl-aux.c
@@ -471,6 +471,33 @@ int QR(KRMAT(x),RMAT(q),RMAT(r)) {
471 OK 471 OK
472} 472}
473 473
474int QRpacked(KRMAT(x),RMAT(qr),RVEC(tau)) {
475 //REQUIRES(xr==rr && xc==rc && qr==qc && xr==qr,BAD_SIZE);
476 DEBUGMSG("QRpacked");
477 KDMVIEW(x);
478 DMVIEW(qr);
479 DVVIEW(tau);
480 //gsl_matrix * a = gsl_matrix_alloc(xr,xc);
481 //gsl_vector * tau = gsl_vector_alloc(MIN(xr,xc));
482 gsl_matrix_memcpy(M(qr),M(x));
483 int res = gsl_linalg_QR_decomp(M(qr),V(tau));
484 CHECK(res,res);
485 OK
486}
487
488
489int QRunpack(KRMAT(xqr),KRVEC(tau),RMAT(q),RMAT(r)) {
490 //REQUIRES(xr==rr && xc==rc && qr==qc && xr==qr,BAD_SIZE);
491 DEBUGMSG("QRunpack");
492 KDMVIEW(xqr);
493 KDVVIEW(tau);
494 DMVIEW(q);
495 DMVIEW(r);
496 gsl_linalg_QR_unpack(M(xqr),V(tau),M(q),M(r));
497 OK
498}
499
500
474int cholR(KRMAT(x),RMAT(l)) { 501int cholR(KRMAT(x),RMAT(l)) {
475 REQUIRES(xr==xc && lr==xr && lr==lc,BAD_SIZE); 502 REQUIRES(xr==xc && lr==xr && lr==lc,BAD_SIZE);
476 DEBUGMSG("cholR"); 503 DEBUGMSG("cholR");
diff --git a/lib/Numeric/GSL/gsl-aux.h b/lib/Numeric/GSL/gsl-aux.h
index 3ccac25..cf22024 100644
--- a/lib/Numeric/GSL/gsl-aux.h
+++ b/lib/Numeric/GSL/gsl-aux.h
@@ -38,6 +38,9 @@ int eigensystemC(KCMAT(x),RVEC(l),CMAT(v));
38 38
39int QR(KRMAT(x),RMAT(q),RMAT(r)); 39int QR(KRMAT(x),RMAT(q),RMAT(r));
40 40
41int QRpacked(KRMAT(x),RMAT(qr),RVEC(tau));
42int QRunpack(KRMAT(qr),KRVEC(tau),RMAT(q),RMAT(r));
43
41int cholR(KRMAT(x),RMAT(l)); 44int cholR(KRMAT(x),RMAT(l));
42 45
43int cholC(KCMAT(x),CMAT(l)); 46int cholC(KCMAT(x),CMAT(l));
diff --git a/lib/Numeric/LinearAlgebra/Algorithms.hs b/lib/Numeric/LinearAlgebra/Algorithms.hs
index 3513b18..5492e07 100644
--- a/lib/Numeric/LinearAlgebra/Algorithms.hs
+++ b/lib/Numeric/LinearAlgebra/Algorithms.hs
@@ -38,17 +38,18 @@ module Numeric.LinearAlgebra.Algorithms (
38 eps, i, 38 eps, i,
39 ctrans, 39 ctrans,
40 Normed(..), NormType(..), 40 Normed(..), NormType(..),
41 GenMat(linearSolveSVD,lu,eigSH') 41 GenMat(linearSolveSVD,lu,eigSH'), unpackQR
42) where 42) where
43 43
44 44
45import Data.Packed.Internal hiding (fromComplex, toComplex, comp, conj) 45import Data.Packed.Internal hiding (fromComplex, toComplex, comp, conj)
46import Data.Packed 46import Data.Packed
47import Numeric.GSL.Matrix(luR,luC,qr) 47import qualified Numeric.GSL.Matrix as GSL
48import Numeric.GSL.Vector 48import Numeric.GSL.Vector
49import Numeric.LinearAlgebra.LAPACK as LAPACK 49import Numeric.LinearAlgebra.LAPACK as LAPACK
50import Complex 50import Complex
51import Numeric.LinearAlgebra.Linear 51import Numeric.LinearAlgebra.Linear
52import Data.List(foldl1')
52 53
53-- | Auxiliary typeclass used to define generic computations for both real and complex matrices. 54-- | Auxiliary typeclass used to define generic computations for both real and complex matrices.
54class (Linear Matrix t) => GenMat t where 55class (Linear Matrix t) => GenMat t where
@@ -59,28 +60,31 @@ class (Linear Matrix t) => GenMat t where
59 eig :: Matrix t -> (Vector (Complex Double), Matrix (Complex Double)) 60 eig :: Matrix t -> (Vector (Complex Double), Matrix (Complex Double))
60 eigSH' :: Matrix t -> (Vector Double, Matrix t) 61 eigSH' :: Matrix t -> (Vector Double, Matrix t)
61 cholSH :: Matrix t -> Matrix t 62 cholSH :: Matrix t -> Matrix t
63 qr :: Matrix t -> (Matrix t, Matrix t)
62 -- | conjugate transpose 64 -- | conjugate transpose
63 ctrans :: Matrix t -> Matrix t 65 ctrans :: Matrix t -> Matrix t
64 66
65instance GenMat Double where 67instance GenMat Double where
66 svd = svdR 68 svd = svdR
67 lu = luR 69 lu = GSL.luR
68 linearSolve = linearSolveR 70 linearSolve = linearSolveR
69 linearSolveSVD = linearSolveSVDR Nothing 71 linearSolveSVD = linearSolveSVDR Nothing
70 ctrans = trans 72 ctrans = trans
71 eig = eigR 73 eig = eigR
72 eigSH' = eigS 74 eigSH' = eigS
73 cholSH = cholS 75 cholSH = cholS
76 qr = GSL.unpackQR . qrR
74 77
75instance GenMat (Complex Double) where 78instance GenMat (Complex Double) where
76 svd = svdC 79 svd = svdC
77 lu = luC 80 lu = GSL.luC
78 linearSolve = linearSolveC 81 linearSolve = linearSolveC
79 linearSolveSVD = linearSolveSVDC Nothing 82 linearSolveSVD = linearSolveSVDC Nothing
80 ctrans = conjTrans 83 ctrans = conjTrans
81 eig = eigC 84 eig = eigC
82 eigSH' = eigH 85 eigSH' = eigH
83 cholSH = cholH 86 cholSH = cholH
87 qr = unpackQR . qrC
84 88
85-- | eigensystem of complex hermitian or real symmetric matrix 89-- | eigensystem of complex hermitian or real symmetric matrix
86eigSH :: GenMat t => Matrix t -> (Vector Double, Matrix t) 90eigSH :: GenMat t => Matrix t -> (Vector Double, Matrix t)
@@ -92,9 +96,6 @@ chol :: GenMat t => Matrix t -> Matrix t
92chol m | m `equal` ctrans m = cholSH m 96chol m | m `equal` ctrans m = cholSH m
93 | otherwise = error "chol requires positive definite complex hermitian or real symmetric matrix" 97 | otherwise = error "chol requires positive definite complex hermitian or real symmetric matrix"
94 98
95qr :: Matrix Double -> (Matrix Double, Matrix Double)
96qr = Numeric.GSL.Matrix.qr
97
98square m = rows m == cols m 99square m = rows m == cols m
99 100
100det :: GenMat t => Matrix t -> t 101det :: GenMat t => Matrix t -> t
@@ -257,3 +258,26 @@ pinvTol t m = v' `mXm` diag s' `mXm` trans u' where
257 d = dim s 258 d = dim s
258 u' = takeColumns d u 259 u' = takeColumns d u
259 v' = takeColumns d v 260 v' = takeColumns d v
261
262---------------------------------------------------------------------
263
264-- many thanks, quickcheck!
265
266haussholder tau v = ident (dim v) `sub` (tau `scale` (w `mXm` ctrans w))
267 where w = asColumn v
268
269unpackQR (pq, tau) = (q,r)
270 where cs = toColumns pq
271 m = rows pq
272 n = cols pq
273 mn = min m n
274 r = fromColumns $ zipWith zt ([m-1, m-2 .. 1] ++ repeat 0) cs
275 vs = zipWith zh [1..mn] cs
276 hs = zipWith haussholder (toList tau) vs
277 q = foldl1' mXm hs
278
279 zh k v = fromList $ replicate (k-1) 0 ++ (1:drop k xs)
280 where xs = toList v
281
282 zt 0 v = v
283 zt k v = join [subVector 0 (dim v - k) v, constant 0 k]
diff --git a/lib/Numeric/LinearAlgebra/LAPACK/lapack-aux.c b/lib/Numeric/LinearAlgebra/LAPACK/lapack-aux.c
index 4ef9a6e..bf0cfba 100644
--- a/lib/Numeric/LinearAlgebra/LAPACK/lapack-aux.c
+++ b/lib/Numeric/LinearAlgebra/LAPACK/lapack-aux.c
@@ -621,7 +621,6 @@ int chol_l_S(KDMAT(a),DMAT(l)) {
621} 621}
622 622
623//////////////////// QR factorization ///////////////////////// 623//////////////////// QR factorization /////////////////////////
624// TO DO: unpack
625 624
626int qr_l_R(KDMAT(a), DVEC(tau), DMAT(r)) { 625int qr_l_R(KDMAT(a), DVEC(tau), DMAT(r)) {
627 integer m = ar; 626 integer m = ar;
@@ -629,7 +628,7 @@ int qr_l_R(KDMAT(a), DVEC(tau), DMAT(r)) {
629 integer mn = MIN(m,n); 628 integer mn = MIN(m,n);
630 REQUIRES(m>=1 && n >=1 && rr== m && rc == n && taun == mn, BAD_SIZE); 629 REQUIRES(m>=1 && n >=1 && rr== m && rc == n && taun == mn, BAD_SIZE);
631 DEBUGMSG("qr_l_R"); 630 DEBUGMSG("qr_l_R");
632 double *WORK = (double*)malloc(m*sizeof(double)); 631 double *WORK = (double*)malloc(n*sizeof(double));
633 CHECK(!WORK,MEM); 632 CHECK(!WORK,MEM);
634 memcpy(rp,ap,m*n*sizeof(double)); 633 memcpy(rp,ap,m*n*sizeof(double));
635 integer res; 634 integer res;
@@ -645,7 +644,7 @@ int qr_l_C(KCMAT(a), CVEC(tau), CMAT(r)) {
645 integer mn = MIN(m,n); 644 integer mn = MIN(m,n);
646 REQUIRES(m>=1 && n >=1 && rr== m && rc == n && taun == mn, BAD_SIZE); 645 REQUIRES(m>=1 && n >=1 && rr== m && rc == n && taun == mn, BAD_SIZE);
647 DEBUGMSG("qr_l_C"); 646 DEBUGMSG("qr_l_C");
648 doublecomplex *WORK = (doublecomplex*)malloc(m*sizeof(doublecomplex)); 647 doublecomplex *WORK = (doublecomplex*)malloc(n*sizeof(doublecomplex));
649 CHECK(!WORK,MEM); 648 CHECK(!WORK,MEM);
650 memcpy(rp,ap,m*n*sizeof(doublecomplex)); 649 memcpy(rp,ap,m*n*sizeof(doublecomplex));
651 integer res; 650 integer res;