1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
|
{-# OPTIONS_GHC -fglasgow-exts #-}
-----------------------------------------------------------------------------
{- |
Module : Numeric.LinearAlgebra.Algorithms
Copyright : (c) Alberto Ruiz 2006-7
License : GPL-style
Maintainer : Alberto Ruiz (aruiz at um dot es)
Stability : provisional
Portability : uses ffi
A generic interface for some common functions. Using it we can write higher level algorithms
and testing properties both for real and complex matrices.
In any case, the specific functions for particular base types can also be explicitly
imported from "Numeric.LinearAlgebra.LAPACK".
-}
-----------------------------------------------------------------------------
module Numeric.LinearAlgebra.Algorithms (
-- * Linear Systems
linearSolve,
inv, pinv,
pinvTol, det, rank, rcond,
-- * Matrix factorizations
-- ** Singular value decomposition
svd,
full, economy, --thin,
-- ** Eigensystems
eig, eigSH,
-- ** QR
qr,
-- ** Cholesky
chol,
-- ** Hessenberg
hess,
-- ** Schur
schur,
-- ** LU
lu,
-- * Matrix functions
expm,
sqrtm,
matFunc,
-- * Nullspace
nullspacePrec,
nullVector,
-- * Norms
Normed(..), NormType(..),
-- * Misc
ctrans,
eps, i,
-- * Util
haussholder,
unpackQR, unpackHess,
Field(linearSolveSVD,eigSH',cholSH)
) where
import Data.Packed.Internal hiding (fromComplex, toComplex, comp, conj, (//))
import Data.Packed
import qualified Numeric.GSL.Matrix as GSL
import Numeric.GSL.Vector
import Numeric.LinearAlgebra.LAPACK as LAPACK
import Complex
import Numeric.LinearAlgebra.Linear
import Data.List(foldl1')
import Data.Array
-- | Auxiliary typeclass used to define generic computations for both real and complex matrices.
class (Normed (Matrix t), Linear Vector t, Linear Matrix t) => Field t where
-- | Singular value decomposition using lapack's dgesvd or zgesvd.
svd :: Matrix t -> (Matrix t, Vector Double, Matrix t)
luPacked :: Matrix t -> (Matrix t, [Int])
-- | Solution of a general linear system (for several right-hand sides) using lapacks' dgesv and zgesv.
-- See also other versions of linearSolve in "Numeric.LinearAlgebra.LAPACK".
linearSolve :: Matrix t -> Matrix t -> Matrix t
linearSolveSVD :: Matrix t -> Matrix t -> Matrix t
-- | Eigenvalues and eigenvectors of a general square matrix using lapack's dgeev or zgeev.
--
-- If @(s,v) = eig m@ then @m \<> v == v \<> diag s@
eig :: Matrix t -> (Vector (Complex Double), Matrix (Complex Double))
-- | Similar to eigSH without checking that the input matrix is hermitian or symmetric.
eigSH' :: Matrix t -> (Vector Double, Matrix t)
-- | Similar to chol without checking that the input matrix is hermitian or symmetric.
cholSH :: Matrix t -> Matrix t
-- | QR factorization using lapack's dgeqr2 or zgeqr2.
--
-- If @(q,r) = qr m@ then @m == q \<> r@, where q is unitary and r is upper triangular.
qr :: Matrix t -> (Matrix t, Matrix t)
-- | Hessenberg factorization using lapack's dgehrd or zgehrd.
--
-- If @(p,h) = hess m@ then @m == p \<> h \<> ctrans p@, where p is unitary
-- and h is in upper Hessenberg form.
hess :: Matrix t -> (Matrix t, Matrix t)
-- | Schur factorization using lapack's dgees or zgees.
--
-- If @(u,s) = schur m@ then @m == u \<> s \<> ctrans u@, where u is unitary
-- and s is a Shur matrix. A complex Schur matrix is upper triangular. A real Schur matrix is
-- upper triangular in 2x2 blocks.
--
-- \"Anything that the Jordan decomposition can do, the Schur decomposition
-- can do better!\" (Van Loan)
schur :: Matrix t -> (Matrix t, Matrix t)
-- | Conjugate transpose.
ctrans :: Matrix t -> Matrix t
instance Field Double where
svd = svdR
luPacked = luR
linearSolve = linearSolveR
linearSolveSVD = linearSolveSVDR Nothing
ctrans = trans
eig = eigR
eigSH' = eigS
cholSH = cholS
qr = GSL.unpackQR . qrR
hess = unpackHess hessR
schur = schurR
instance Field (Complex Double) where
svd = svdC
luPacked = luC
linearSolve = linearSolveC
linearSolveSVD = linearSolveSVDC Nothing
ctrans = conj . trans
eig = eigC
eigSH' = eigH
cholSH = cholH
qr = unpackQR . qrC
hess = unpackHess hessC
schur = schurC
-- | Eigenvalues and Eigenvectors of a complex hermitian or real symmetric matrix using lapack's dsyev or zheev.
--
-- If @(s,v) = eigSH m@ then @m == v \<> diag s \<> ctrans v@
eigSH :: Field t => Matrix t -> (Vector Double, Matrix t)
eigSH m | m `equal` ctrans m = eigSH' m
| otherwise = error "eigSH requires complex hermitian or real symmetric matrix"
-- | Cholesky factorization of a positive definite hermitian or symmetric matrix using lapack's dpotrf or zportrf.
--
-- If @c = chol m@ then @m == c \<> ctrans c@.
chol :: Field t => Matrix t -> Matrix t
chol m | m `equal` ctrans m = cholSH m
| otherwise = error "chol requires positive definite complex hermitian or real symmetric matrix"
square m = rows m == cols m
-- | determinant of a square matrix, computed from the LU decomposition.
det :: Field t => Matrix t -> t
det m | square m = s * (product $ toList $ takeDiag $ lu)
| otherwise = error "det of nonsquare matrix"
where (lu,perm) = luPacked m
s = signlp (rows m) perm
-- | LU factorization of a general matrix using lapack's dgetrf or zgetrf.
--
-- If @(l,u,p,s) = lu m@ then @m == p \<> l \<> u@, where l is lower triangular,
-- u is upper triangular, p is a permutation matrix and s is the signature of the permutation.
lu :: Field t => Matrix t -> (Matrix t, Matrix t, Matrix t, t)
lu = luFact . luPacked
-- | Inverse of a square matrix using lapacks' dgesv and zgesv.
inv :: Field t => Matrix t -> Matrix t
inv m | square m = m `linearSolve` ident (rows m)
| otherwise = error "inv of nonsquare matrix"
-- | Pseudoinverse of a general matrix using lapack's dgelss or zgelss.
pinv :: Field t => Matrix t -> Matrix t
pinv m = linearSolveSVD m (ident (rows m))
-- | A version of 'svd' which returns an appropriate diagonal matrix with the singular values.
--
-- If @(u,d,v) = full svd m@ then @m == u \<> d \<> trans v@.
full :: Element t
=> (Matrix t -> (Matrix t, Vector Double, Matrix t)) -> Matrix t -> (Matrix t, Matrix Double, Matrix t)
full svd' m = (u, d ,v) where
(u,s,v) = svd' m
d = diagRect s r c
r = rows m
c = cols m
-- | A version of 'svd' which returns only the nonzero singular values and the corresponding rows and columns of the rotations.
--
-- If @(u,s,v) = economy svd m@ then @m == u \<> diag s \<> trans v@.
economy :: Element t
=> (Matrix t -> (Matrix t, Vector Double, Matrix t)) -> Matrix t -> (Matrix t, Vector Double, Matrix t)
economy svd' m = (u', subVector 0 d s, v') where
(u,s,v) = svd' m
sl@(g:_) = toList s
s' = fromList . filter (>tol) $ sl
t = 1
tol = (fromIntegral (max r c) * g * t * eps)
r = rows m
c = cols m
d = dim s'
u' = takeColumns d u
v' = takeColumns d v
-- | The machine precision of a Double: @eps = 2.22044604925031e-16@ (the value used by GNU-Octave).
eps :: Double
eps = 2.22044604925031e-16
-- | The imaginary unit: @i = 0.0 :+ 1.0@
i :: Complex Double
i = 0:+1
-- matrix product
mXm :: (Num t, Field t) => Matrix t -> Matrix t -> Matrix t
mXm = multiply
-- matrix - vector product
mXv :: (Num t, Field t) => Matrix t -> Vector t -> Vector t
mXv m v = flatten $ m `mXm` (asColumn v)
-- vector - matrix product
vXm :: (Num t, Field t) => Vector t -> Matrix t -> Vector t
vXm v m = flatten $ (asRow v) `mXm` m
---------------------------------------------------------------------------
norm2 :: Vector Double -> Double
norm2 = toScalarR Norm2
norm1 :: Vector Double -> Double
norm1 = toScalarR AbsSum
data NormType = Infinity | PNorm1 | PNorm2 -- PNorm Int
pnormRV PNorm2 = norm2
pnormRV PNorm1 = norm1
pnormRV Infinity = vectorMax . vectorMapR Abs
--pnormRV _ = error "pnormRV not yet defined"
pnormCV PNorm2 = norm2 . asReal
pnormCV PNorm1 = norm1 . liftVector magnitude
pnormCV Infinity = vectorMax . liftVector magnitude
--pnormCV _ = error "pnormCV not yet defined"
pnormRM PNorm2 m = head (toList s) where (_,s,_) = svdR m
pnormRM PNorm1 m = vectorMax $ constant 1 (rows m) `vXm` liftMatrix (vectorMapR Abs) m
pnormRM Infinity m = vectorMax $ liftMatrix (vectorMapR Abs) m `mXv` constant 1 (cols m)
--pnormRM _ _ = error "p norm not yet defined"
pnormCM PNorm2 m = head (toList s) where (_,s,_) = svdC m
pnormCM PNorm1 m = vectorMax $ constant 1 (rows m) `vXm` liftMatrix (liftVector magnitude) m
pnormCM Infinity m = vectorMax $ liftMatrix (liftVector magnitude) m `mXv` constant 1 (cols m)
--pnormCM _ _ = error "p norm not yet defined"
-- | Objects which have a p-norm.
-- Using it you can define convenient shortcuts:
--
-- @norm2 x = pnorm PNorm2 x@
--
-- @frobenius m = norm2 . flatten $ m@
class Normed t where
pnorm :: NormType -> t -> Double
instance Normed (Vector Double) where
pnorm = pnormRV
instance Normed (Vector (Complex Double)) where
pnorm = pnormCV
instance Normed (Matrix Double) where
pnorm = pnormRM
instance Normed (Matrix (Complex Double)) where
pnorm = pnormCM
-----------------------------------------------------------------------
-- | The nullspace of a matrix from its SVD decomposition.
nullspacePrec :: Field t
=> Double -- ^ relative tolerance in 'eps' units
-> Matrix t -- ^ input matrix
-> [Vector t] -- ^ list of unitary vectors spanning the nullspace
nullspacePrec t m = ns where
(_,s,v) = svd m
sl@(g:_) = toList s
tol = (fromIntegral (max (rows m) (cols m)) * g * t * eps)
rank' = length (filter (> g*tol) sl)
ns = drop rank' $ toRows $ ctrans v
-- | The nullspace of a matrix, assumed to be one-dimensional, with default tolerance (shortcut for @last . nullspacePrec 1@).
nullVector :: Field t => Matrix t -> Vector t
nullVector = last . nullspacePrec 1
------------------------------------------------------------------------
{- Pseudoinverse of a real matrix with the desired tolerance, expressed as a
multiplicative factor of the default tolerance used by GNU-Octave (see 'pinv').
@\> let m = 'fromLists' [[1,0, 0]
,[0,1, 0]
,[0,0,1e-10]]
\
\> 'pinv' m
1. 0. 0.
0. 1. 0.
0. 0. 10000000000.
\
\> pinvTol 1E8 m
1. 0. 0.
0. 1. 0.
0. 0. 1.@
-}
--pinvTol :: Double -> Matrix Double -> Matrix Double
pinvTol t m = v' `mXm` diag s' `mXm` trans u' where
(u,s,v) = svdR m
sl@(g:_) = toList s
s' = fromList . map rec $ sl
rec x = if x < g*tol then 1 else 1/x
tol = (fromIntegral (max r c) * g * t * eps)
r = rows m
c = cols m
d = dim s
u' = takeColumns d u
v' = takeColumns d v
---------------------------------------------------------------------
-- many thanks, quickcheck!
haussholder :: (Field a) => a -> Vector a -> Matrix a
haussholder tau v = ident (dim v) `sub` (tau `scale` (w `mXm` ctrans w))
where w = asColumn v
zh k v = fromList $ replicate (k-1) 0 ++ (1:drop k xs)
where xs = toList v
zt 0 v = v
zt k v = join [subVector 0 (dim v - k) v, constant 0 k]
unpackQR :: (Field t) => (Matrix t, Vector t) -> (Matrix t, Matrix t)
unpackQR (pq, tau) = (q,r)
where cs = toColumns pq
m = rows pq
n = cols pq
mn = min m n
r = fromColumns $ zipWith zt ([m-1, m-2 .. 1] ++ repeat 0) cs
vs = zipWith zh [1..mn] cs
hs = zipWith haussholder (toList tau) vs
q = foldl1' mXm hs
unpackHess :: (Field t) => (Matrix t -> (Matrix t,Vector t)) -> Matrix t -> (Matrix t, Matrix t)
unpackHess hf m
| rows m == 1 = ((1><1)[1],m)
| otherwise = (uH . hf) m
uH (pq, tau) = (p,h)
where cs = toColumns pq
m = rows pq
n = cols pq
mn = min m n
h = fromColumns $ zipWith zt ([m-2, m-3 .. 1] ++ repeat 0) cs
vs = zipWith zh [2..mn] cs
hs = zipWith haussholder (toList tau) vs
p = foldl1' mXm hs
--------------------------------------------------------------------------
-- | Reciprocal of the 2-norm condition number of a matrix, computed from the SVD.
rcond :: Field t => Matrix t -> Double
rcond m = last s / head s
where (_,s',_) = svd m
s = toList s'
-- | Number of linearly independent rows or columns.
rank :: Field t => Matrix t -> Int
rank m | pnorm PNorm1 m < eps = 0
| otherwise = dim s where (_,s,_) = economy svd m
{-
expm' m = case diagonalize (complex m) of
Just (l,v) -> v `mXm` diag (exp l) `mXm` inv v
Nothing -> error "Sorry, expm not yet implemented for non-diagonalizable matrices"
where exp = vectorMapC Exp
-}
diagonalize m = if rank v == n
then Just (l,v)
else Nothing
where n = rows m
(l,v) = if m `equal` ctrans m
then let (l',v') = eigSH m in (real l', v')
else eig m
-- | Generic matrix functions for diagonalizable matrices. For instance:
--
-- @logm = matFunc log@
--
matFunc :: Field t => (Complex Double -> Complex Double) -> Matrix t -> Matrix (Complex Double)
matFunc f m = case diagonalize (complex m) of
Just (l,v) -> v `mXm` diag (liftVector f l) `mXm` inv v
Nothing -> error "Sorry, matFunc requires a diagonalizable matrix"
--------------------------------------------------------------
golubeps :: Integer -> Integer -> Double
golubeps p q = a * fromIntegral b / fromIntegral c where
a = 2^^(3-p-q)
b = fact p * fact q
c = fact (p+q) * fact (p+q+1)
fact n = product [1..n]
epslist = [ (fromIntegral k, golubeps k k) | k <- [1..]]
geps delta = head [ k | (k,g) <- epslist, g<delta]
expGolub m = iterate msq f !! j
where j = max 0 $ floor $ log2 $ pnorm Infinity m
log2 x = log x / log 2
a = m */ fromIntegral ((2::Int)^j)
q = geps eps -- 7 steps
eye = ident (rows m)
work (k,c,x,n,d) = (k',c',x',n',d')
where k' = k+1
c' = c * fromIntegral (q-k+1) / fromIntegral ((2*q-k+1)*k)
x' = a <> x
n' = n |+| (c' .* x')
d' = d |+| (((-1)^k * c') .* x')
(_,_,_,nf,df) = iterate work (1,1,eye,eye,eye) !! q
f = linearSolve df nf
msq x = x <> x
(<>) = multiply
v */ x = scale (recip x) v
(.*) = scale
(|+|) = add
{- | Matrix exponential. It uses a direct translation of Algorithm 11.3.1 in Golub & Van Loan,
based on a scaled Pade approximation.
-}
expm :: Field t => Matrix t -> Matrix t
expm = expGolub
--------------------------------------------------------------
{- | Matrix square root. Currently it uses a simple iterative algorithm described in Wikipedia.
It only works with invertible matrices that have a real solution. For diagonalizable matrices you can try @matFunc sqrt@.
@m = (2><2) [4,9
,0,4] :: Matrix Double@
@\>sqrtm m
(2><2)
[ 2.0, 2.25
, 0.0, 2.0 ]@
-}
sqrtm :: Field t => Matrix t -> Matrix t
sqrtm = sqrtmInv
sqrtmInv x = fst $ fixedPoint $ iterate f (x, ident (rows x))
where fixedPoint (a:b:rest) | pnorm PNorm1 (fst a |-| fst b) < eps = a
| otherwise = fixedPoint (b:rest)
fixedPoint _ = error "fixedpoint with impossible inputs"
f (y,z) = (0.5 .* (y |+| inv z),
0.5 .* (inv y |+| z))
(.*) = scale
(|+|) = add
(|-|) = sub
------------------------------------------------------------------
signlp r vals = foldl f 1 (zip [0..r-1] vals)
where f s (a,b) | a /= b = -s
| otherwise = s
swap (arr,s) (a,b) | a /= b = (arr // [(a, arr!b),(b,arr!a)],-s)
| otherwise = (arr,s)
fixPerm r vals = (fromColumns $ elems res, sign)
where v = [0..r-1]
s = toColumns (ident r)
(res,sign) = foldl swap (listArray (0,r-1) s, 1) (zip v vals)
triang r c h v = reshape c $ fromList [el i j | i<-[0..r-1], j<-[0..c-1]]
where el i j = if j-i>=h then v else 1 - v
luFact (lu,perm) | r <= c = (l ,u ,p, s)
| otherwise = (l',u',p, s)
where
r = rows lu
c = cols lu
tu = triang r c 0 1
tl = triang r c 0 0
l = takeColumns r (lu |*| tl) |+| diagRect (constant 1 r) r r
u = lu |*| tu
(p,s) = fixPerm r perm
l' = (lu |*| tl) |+| diagRect (constant 1 c) r c
u' = takeRows c (lu |*| tu)
(|+|) = add
(|*|) = mul
|