{-# LANGUAGE FlexibleContexts, FlexibleInstances #-}
{-# LANGUAGE CPP #-}
{-# LANGUAGE MultiParamTypeClasses #-}
{-# LANGUAGE UndecidableInstances #-}
{-# LANGUAGE TypeFamilies #-}
{-# OPTIONS_GHC -fno-warn-missing-signatures #-}
module Internal.Algorithms (
module Internal.Algorithms,
UpLo(..)
) where
import Internal.Vector
import Internal.Matrix
import Internal.Element
import Internal.Conversion
import Internal.LAPACK
import Internal.Numeric
import Data.List(foldl1')
import qualified Data.Array as A
import qualified Data.Vector.Storable as Vector
import Internal.ST
import Internal.Vectorized(range)
import Control.DeepSeq
class (Numeric t,
Convert t,
Normed Matrix t,
Normed Vector t,
Floating t,
Linear t Vector,
Linear t Matrix,
Additive (Vector t),
Additive (Matrix t),
RealOf t ~ Double) => Field t where
svd' :: Matrix t -> (Matrix t, Vector Double, Matrix t)
thinSVD' :: Matrix t -> (Matrix t, Vector Double, Matrix t)
sv' :: Matrix t -> Vector Double
luPacked' :: Matrix t -> (Matrix t, [Int])
luSolve' :: (Matrix t, [Int]) -> Matrix t -> Matrix t
mbLinearSolve' :: Matrix t -> Matrix t -> Maybe (Matrix t)
linearSolve' :: Matrix t -> Matrix t -> Matrix t
cholSolve' :: Matrix t -> Matrix t -> Matrix t
triSolve' :: UpLo -> Matrix t -> Matrix t -> Matrix t
triDiagSolve' :: Vector t -> Vector t -> Vector t -> Matrix t -> Matrix t
ldlPacked' :: Matrix t -> (Matrix t, [Int])
ldlSolve' :: (Matrix t, [Int]) -> Matrix t -> Matrix t
linearSolveSVD' :: Matrix t -> Matrix t -> Matrix t
linearSolveLS' :: Matrix t -> Matrix t -> Matrix t
eig' :: Matrix t -> (Vector (Complex Double), Matrix (Complex Double))
eigSH'' :: Matrix t -> (Vector Double, Matrix t)
eigOnly :: Matrix t -> Vector (Complex Double)
eigOnlySH :: Matrix t -> Vector Double
cholSH' :: Matrix t -> Matrix t
mbCholSH' :: Matrix t -> Maybe (Matrix t)
qr' :: Matrix t -> (Matrix t, Vector t)
qrgr' :: Int -> (Matrix t, Vector t) -> Matrix t
hess' :: Matrix t -> (Matrix t, Matrix t)
schur' :: Matrix t -> (Matrix t, Matrix t)
instance Field Double where
svd' = svdRd
thinSVD' = thinSVDRd
sv' = svR
luPacked' = luR
luSolve' (l_u,perm) = lusR l_u perm
linearSolve' = linearSolveR
mbLinearSolve' = mbLinearSolveR
cholSolve' = cholSolveR
triSolve' = triSolveR
triDiagSolve' = triDiagSolveR
linearSolveLS' = linearSolveLSR
linearSolveSVD' = linearSolveSVDR Nothing
eig' = eigR
eigSH'' = eigS
eigOnly = eigOnlyR
eigOnlySH = eigOnlyS
cholSH' = cholS
mbCholSH' = mbCholS
qr' = qrR
qrgr' = qrgrR
hess' = unpackHess hessR
schur' = schurR
ldlPacked' = ldlR
ldlSolve'= uncurry ldlsR
instance Field (Complex Double) where
#ifdef NOZGESDD
svd' = svdC
thinSVD' = thinSVDC
#else
svd' = svdCd
thinSVD' = thinSVDCd
#endif
sv' = svC
luPacked' = luC
luSolve' (l_u,perm) = lusC l_u perm
linearSolve' = linearSolveC
mbLinearSolve' = mbLinearSolveC
cholSolve' = cholSolveC
triSolve' = triSolveC
triDiagSolve' = triDiagSolveC
linearSolveLS' = linearSolveLSC
linearSolveSVD' = linearSolveSVDC Nothing
eig' = eigC
eigOnly = eigOnlyC
eigSH'' = eigH
eigOnlySH = eigOnlyH
cholSH' = cholH
mbCholSH' = mbCholH
qr' = qrC
qrgr' = qrgrC
hess' = unpackHess hessC
schur' = schurC
ldlPacked' = ldlC
ldlSolve' = uncurry ldlsC
square m = rows m == cols m
vertical m = rows m >= cols m
exactHermitian m = m `equal` ctrans m
svd :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)
svd = {-# SCC "svd" #-} g . svd'
where
g (u,s,v) = (u,s,tr v)
thinSVD :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)
thinSVD = {-# SCC "thinSVD" #-} g . thinSVD'
where
g (u,s,v) = (u,s,tr v)
singularValues :: Field t => Matrix t -> Vector Double
singularValues = {-# SCC "singularValues" #-} sv'
fullSVD :: Field t => Matrix t -> (Matrix t, Matrix Double, Matrix t)
fullSVD m = (u,d,v) where
(u,s,v) = svd m
d = diagRect 0 s r c
r = rows m
c = cols m
compactSVD :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)
compactSVD = compactSVDTol 1
compactSVDTol :: Field t => Double -> Matrix t -> (Matrix t, Vector Double, Matrix t)
compactSVDTol r m = (u', subVector 0 d s, v') where
(u,s,v) = thinSVD m
d = rankSVD (r*eps) m s `max` 1
u' = takeColumns d u
v' = takeColumns d v
rightSV :: Field t => Matrix t -> (Vector Double, Matrix t)
rightSV m | vertical m = let (_,s,v) = thinSVD m in (s,v)
| otherwise = let (_,s,v) = svd m in (s,v)
leftSV :: Field t => Matrix t -> (Matrix t, Vector Double)
leftSV m | vertical m = let (u,s,_) = svd m in (u,s)
| otherwise = let (u,s,_) = thinSVD m in (u,s)
data LU t = LU (Matrix t) [Int] deriving Show
instance (NFData t, Numeric t) => NFData (LU t)
where
rnf (LU m _) = rnf m
luPacked :: Field t => Matrix t -> LU t
luPacked x = {-# SCC "luPacked" #-} LU m p
where
(m,p) = luPacked' x
luSolve :: Field t => LU t -> Matrix t -> Matrix t
luSolve (LU m p) = {-# SCC "luSolve" #-} luSolve' (m,p)
linearSolve :: Field t => Matrix t -> Matrix t -> Matrix t
linearSolve = {-# SCC "linearSolve" #-} linearSolve'
mbLinearSolve :: Field t => Matrix t -> Matrix t -> Maybe (Matrix t)
mbLinearSolve = {-# SCC "linearSolve" #-} mbLinearSolve'
cholSolve
:: Field t
=> Matrix t
-> Matrix t
-> Matrix t
cholSolve = {-# SCC "cholSolve" #-} cholSolve'
triSolve
:: Field t
=> UpLo
-> Matrix t
-> Matrix t
-> Matrix t
triSolve = {-# SCC "triSolve" #-} triSolve'
triDiagSolve
:: Field t
=> Vector t
-> Vector t
-> Vector t
-> Matrix t
-> Matrix t
triDiagSolve = {-# SCC "triDiagSolve" #-} triDiagSolve'
linearSolveSVD :: Field t => Matrix t -> Matrix t -> Matrix t
linearSolveSVD = {-# SCC "linearSolveSVD" #-} linearSolveSVD'
linearSolveLS :: Field t => Matrix t -> Matrix t -> Matrix t
linearSolveLS = {-# SCC "linearSolveLS" #-} linearSolveLS'
data LDL t = LDL (Matrix t) [Int] deriving Show
instance (NFData t, Numeric t) => NFData (LDL t)
where
rnf (LDL m _) = rnf m
ldlPackedSH :: Field t => Matrix t -> LDL t
ldlPackedSH x = {-# SCC "ldlPacked" #-} LDL m p
where
(m,p) = ldlPacked' x
ldlPacked :: Field t => Herm t -> LDL t
ldlPacked (Herm m) = ldlPackedSH m
ldlSolve :: Field t => LDL t -> Matrix t -> Matrix t
ldlSolve (LDL m p) = {-# SCC "ldlSolve" #-} ldlSolve' (m,p)
eig :: Field t => Matrix t -> (Vector (Complex Double), Matrix (Complex Double))
eig = {-# SCC "eig" #-} eig'
eigenvalues :: Field t => Matrix t -> Vector (Complex Double)
eigenvalues = {-# SCC "eigenvalues" #-} eigOnly
eigSH' :: Field t => Matrix t -> (Vector Double, Matrix t)
eigSH' = {-# SCC "eigSH'" #-} eigSH''
eigenvaluesSH' :: Field t => Matrix t -> Vector Double
eigenvaluesSH' = {-# SCC "eigenvaluesSH'" #-} eigOnlySH
eigSH :: Field t => Herm t -> (Vector Double, Matrix t)
eigSH (Herm m) = eigSH' m
eigenvaluesSH :: Field t => Herm t -> Vector Double
eigenvaluesSH (Herm m) = eigenvaluesSH' m
data QR t = QR (Matrix t) (Vector t)
instance (NFData t, Numeric t) => NFData (QR t)
where
rnf (QR m _) = rnf m
qr :: Field t => Matrix t -> (Matrix t, Matrix t)
qr = {-# SCC "qr" #-} unpackQR . qr'
thinQR :: Field t => Matrix t -> (Matrix t, Matrix t)
thinQR = {-# SCC "thinQR" #-} thinUnpackQR . qr'
qrRaw :: Field t => Matrix t -> QR t
qrRaw m = QR x v
where
(x,v) = qr' m
qrgr :: Field t => Int -> QR t -> Matrix t
qrgr n (QR a t)
| dim t > min (cols a) (rows a) || n < 0 || n > dim t = error "qrgr expects k <= min(rows,cols)"
| otherwise = qrgr' n (a,t)
rq :: Field t => Matrix t -> (Matrix t, Matrix t)
rq = {-# SCC "rq" #-} rqFromQR qr
thinRQ :: Field t => Matrix t -> (Matrix t, Matrix t)
thinRQ = {-# SCC "thinQR" #-} rqFromQR thinQR
rqFromQR :: Field t => (Matrix t -> (Matrix t, Matrix t)) -> Matrix t -> (Matrix t, Matrix t)
rqFromQR qr0 m = (r,q) where
(q',r') = qr0 $ trans $ rev1 m
r = rev2 (trans r')
q = rev2 (trans q')
rev1 = flipud . fliprl
rev2 = fliprl . flipud
hess :: Field t => Matrix t -> (Matrix t, Matrix t)
hess = hess'
schur :: Field t => Matrix t -> (Matrix t, Matrix t)
schur = schur'
mbCholSH :: Field t => Matrix t -> Maybe (Matrix t)
mbCholSH = {-# SCC "mbCholSH" #-} mbCholSH'
cholSH :: Field t => Matrix t -> Matrix t
cholSH = cholSH'
chol :: Field t => Herm t -> Matrix t
chol (Herm m) = {-# SCC "chol" #-} cholSH' m
mbChol :: Field t => Herm t -> Maybe (Matrix t)
mbChol (Herm m) = {-# SCC "mbChol" #-} mbCholSH' m
invlndet :: Field t
=> Matrix t
-> (Matrix t, (t, t))
invlndet m | square m = (im,(ladm,sdm))
| otherwise = error $ "invlndet of nonsquare "++ shSize m ++ " matrix"
where
lp@(LU lup perm) = luPacked m
s = signlp (rows m) perm
dg = toList $ takeDiag $ lup
ladm = sum $ map (log.abs) dg
sdm = s* product (map signum dg)
im = luSolve lp (ident (rows m))
det :: Field t => Matrix t -> t
det m | square m = {-# SCC "det" #-} s * (product $ toList $ takeDiag $ lup)
| otherwise = error $ "det of nonsquare "++ shSize m ++ " matrix"
where
LU lup perm = luPacked m
s = signlp (rows m) perm
lu :: Field t => Matrix t -> (Matrix t, Matrix t, Matrix t, t)
lu = luFact . luPacked
inv :: Field t => Matrix t -> Matrix t
inv m | square m = m `linearSolve` ident (rows m)
| otherwise = error $ "inv of nonsquare "++ shSize m ++ " matrix"
pinv :: Field t => Matrix t -> Matrix t
pinv = pinvTol 1
pinvTol :: Field t => Double -> Matrix t -> Matrix t
pinvTol t m = v' `mXm` diag s' `mXm` ctrans u' where
(u,s,v) = thinSVD m
sl@(g:_) = toList s
s' = real . fromList . map rec $ sl
rec x = if x <= g*tol then x 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
rankSVD :: Element t
=> Double
-> Matrix t
-> Vector Double
-> Int
rankSVD teps m s = ranksv teps (max (rows m) (cols m)) (toList s)
ranksv :: Double
-> Int
-> [Double]
-> Int
ranksv teps maxdim s = k where
g = maximum s
tol = fromIntegral maxdim * g * teps
s' = filter (>tol) s
k = if g > teps then length s' else 0
eps :: Double
eps = 2.22044604925031e-16
peps :: RealFloat x => x
peps = x where x = 2.0 ** fromIntegral (1 - floatDigits x)
nullspaceSVD :: Field t
=> Either Double Int
-> Matrix t
-> (Vector Double, Matrix t)
-> Matrix t
nullspaceSVD hint a (s,v) = vs where
tol = case hint of
Left t -> t
_ -> eps
k = case hint of
Right t -> t
_ -> rankSVD tol a s
vs = dropColumns k v
nullspacePrec :: Field t
=> Double
-> Matrix t
-> [Vector t]
nullspacePrec t m = toColumns $ nullspaceSVD (Left (t*eps)) m (rightSV m)
nullVector :: Field t => Matrix t -> Vector t
nullVector = last . nullspacePrec 1
orthSVD :: Field t
=> Either Double Int
-> Matrix t
-> (Matrix t, Vector Double)
-> Matrix t
orthSVD hint a (v,s) = vs where
tol = case hint of
Left t -> t
_ -> eps
k = case hint of
Right t -> t
_ -> rankSVD tol a s
vs = takeColumns k v
orth :: Field t => Matrix t -> [Vector t]
orth m = take r $ toColumns u
where
(u,s,_) = compactSVD m
r = ranksv eps (max (rows m) (cols m)) (toList s)
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 = vjoin [subVector 0 (dim v - k) v, konst' 0 k]
unpackQR :: (Field t) => (Matrix t, Vector t) -> (Matrix t, Matrix t)
unpackQR (pq, tau) = {-# SCC "unpackQR" #-} (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
thinUnpackQR :: (Field t) => (Matrix t, Vector t) -> (Matrix t, Matrix t)
thinUnpackQR (pq, tau) = (q, r)
where mn = uncurry min $ size pq
q = qrgr mn $ QR pq tau
r = fromRows $ zipWith (\i v -> Vector.replicate i 0 Vector.++ Vector.drop i v) [0..mn-1] (toRows pq)
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
rcond :: Field t => Matrix t -> Double
rcond m = last s / head s
where s = toList (singularValues m)
rank :: Field t => Matrix t -> Int
rank m = rankSVD eps m (singularValues m)
diagonalize m = if rank v == n
then Just (l,v)
else Nothing
where n = rows m
(l,v) = if exactHermitian m
then let (l',v') = eigSH (trustSym m) in (real l', v')
else eig m
matFunc :: (Complex Double -> Complex Double) -> Matrix (Complex Double) -> Matrix (Complex Double)
matFunc f m = case diagonalize m of
Just (l,v) -> v `mXm` diag (mapVector 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 :: [(Int,Double)]
epslist = [ (fromIntegral k, golubeps k k) | k <- [1..]]
geps delta = head [ k | (k,g) <- epslist, g<delta]
expm :: Field t => Matrix t -> Matrix t
expm = expGolub
expGolub :: Field t => Matrix t -> Matrix t
expGolub m = iterate msq f !! j
where j = max 0 $ floor $ logBase 2 $ pnorm Infinity m
a = m */ fromIntegral ((2::Int)^j)
q = geps eps
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
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) < peps = 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
fixPerm r vals = (fromColumns $ A.elems res, sign)
where
v = [0..r-1]
t = toColumns (ident r)
(res,sign) = foldl swap (A.listArray (0,r-1) t, 1) (zip v vals)
swap (arr,s) (a,b)
| a /= b = (arr A.// [(a, arr A.! b),(b,arr A.! a)],-s)
| otherwise = (arr,s)
fixPerm' :: [Int] -> Vector I
fixPerm' s = res $ mutable f s0
where
s0 = reshape 1 (range (length s))
res = flatten . fst
swap m i j = rowOper (SWAP i j AllCols) m
f :: (Num t, Element t) => (Int, Int) -> STMatrix s t -> ST s ()
f _ p = sequence_ $ zipWith (swap p) [0..] s
triang r c h v = (r><c) [el s t | s<-[0..r-1], t<-[0..c-1]]
where el p q = if q-p>=h then v else 1 - v
luFact :: Numeric t => LU t -> (Matrix t, Matrix t, Matrix t, t)
luFact (LU l_u perm)
| r <= c = (l ,u ,p, s)
| otherwise = (l',u',p, s)
where
r = rows l_u
c = cols l_u
tu = triang r c 0 1
tl = triang r c 0 0
l = takeColumns r (l_u |*| tl) |+| diagRect 0 (konst' 1 r) r r
u = l_u |*| tu
(p,s) = fixPerm r perm
l' = (l_u |*| tl) |+| diagRect 0 (konst' 1 c) r c
u' = takeRows c (l_u |*| tu)
(|+|) = add
(|*|) = mul
data NormType = Infinity | PNorm1 | PNorm2 | Frobenius
class (RealFloat (RealOf t)) => Normed c t where
pnorm :: NormType -> c t -> RealOf t
instance Normed Vector Double where
pnorm PNorm1 = norm1
pnorm PNorm2 = norm2
pnorm Infinity = normInf
pnorm Frobenius = norm2
instance Normed Vector (Complex Double) where
pnorm PNorm1 = norm1
pnorm PNorm2 = norm2
pnorm Infinity = normInf
pnorm Frobenius = pnorm PNorm2
instance Normed Vector Float where
pnorm PNorm1 = norm1
pnorm PNorm2 = norm2
pnorm Infinity = normInf
pnorm Frobenius = pnorm PNorm2
instance Normed Vector (Complex Float) where
pnorm PNorm1 = norm1
pnorm PNorm2 = norm2
pnorm Infinity = normInf
pnorm Frobenius = pnorm PNorm2
instance Normed Matrix Double where
pnorm PNorm1 = maximum . map (pnorm PNorm1) . toColumns
pnorm PNorm2 = (@>0) . singularValues
pnorm Infinity = pnorm PNorm1 . trans
pnorm Frobenius = pnorm PNorm2 . flatten
instance Normed Matrix (Complex Double) where
pnorm PNorm1 = maximum . map (pnorm PNorm1) . toColumns
pnorm PNorm2 = (@>0) . singularValues
pnorm Infinity = pnorm PNorm1 . trans
pnorm Frobenius = pnorm PNorm2 . flatten
instance Normed Matrix Float where
pnorm PNorm1 = maximum . map (pnorm PNorm1) . toColumns
pnorm PNorm2 = realToFrac . (@>0) . singularValues . double
pnorm Infinity = pnorm PNorm1 . trans
pnorm Frobenius = pnorm PNorm2 . flatten
instance Normed Matrix (Complex Float) where
pnorm PNorm1 = maximum . map (pnorm PNorm1) . toColumns
pnorm PNorm2 = realToFrac . (@>0) . singularValues . double
pnorm Infinity = pnorm PNorm1 . trans
pnorm Frobenius = pnorm PNorm2 . flatten
relativeError' :: (Normed c t, Container c t) => c t -> c t -> Int
relativeError' x y = dig (norm (x `sub` y) / norm x)
where norm = pnorm Infinity
dig r = round $ -logBase 10 (realToFrac r :: Double)
relativeError :: Num a => (a -> Double) -> a -> a -> Double
relativeError norm a b = r
where
na = norm a
nb = norm b
nab = norm (a-b)
mx = max na nb
mn = min na nb
r = if mn < peps
then mx
else nab/mx
geigSH :: Field t
=> Herm t
-> Herm t
-> (Vector Double, Matrix t)
geigSH (Herm a) (Herm b) = geigSH' a b
geigSH' :: Field t
=> Matrix t
-> Matrix t
-> (Vector Double, Matrix t)
geigSH' a b = (l,v')
where
u = cholSH b
iu = inv u
c = ctrans iu <> a <> iu
(l,v) = eigSH' c
v' = iu <> v
(<>) = mXm
newtype Herm t = Herm (Matrix t) deriving Show
instance (NFData t, Numeric t) => NFData (Herm t)
where
rnf (Herm m) = rnf m
unSym :: Herm t -> Matrix t
unSym (Herm x) = x
sym :: Field t => Matrix t -> Herm t
sym x = Herm (scale 0.5 (tr x `add` x))
mTm :: Numeric t => Matrix t -> Herm t
mTm x = Herm (tr x `mXm` x)
instance Field t => Linear t Herm where
scale x (Herm m) = Herm (scale x m)
instance Field t => Additive (Herm t) where
add (Herm a) (Herm b) = Herm (a `add` b)
trustSym :: Matrix t -> Herm t
trustSym x = (Herm x)