diff options
Diffstat (limited to 'packages/base/src/Internal')
-rw-r--r-- | packages/base/src/Internal/Algorithms.hs | 89 | ||||
-rw-r--r-- | packages/base/src/Internal/Modular.hs | 2 | ||||
-rw-r--r-- | packages/base/src/Internal/Util.hs | 2 |
3 files changed, 53 insertions, 40 deletions
diff --git a/packages/base/src/Internal/Algorithms.hs b/packages/base/src/Internal/Algorithms.hs index d2f17f4..c4f1a60 100644 --- a/packages/base/src/Internal/Algorithms.hs +++ b/packages/base/src/Internal/Algorithms.hs | |||
@@ -4,12 +4,6 @@ | |||
4 | {-# LANGUAGE UndecidableInstances #-} | 4 | {-# LANGUAGE UndecidableInstances #-} |
5 | {-# LANGUAGE TypeFamilies #-} | 5 | {-# LANGUAGE TypeFamilies #-} |
6 | 6 | ||
7 | {-# LANGUAGE FlexibleContexts, FlexibleInstances #-} | ||
8 | {-# LANGUAGE CPP #-} | ||
9 | {-# LANGUAGE MultiParamTypeClasses #-} | ||
10 | {-# LANGUAGE UndecidableInstances #-} | ||
11 | {-# LANGUAGE TypeFamilies #-} | ||
12 | |||
13 | ----------------------------------------------------------------------------- | 7 | ----------------------------------------------------------------------------- |
14 | {- | | 8 | {- | |
15 | Module : Internal.Algorithms | 9 | Module : Internal.Algorithms |
@@ -376,8 +370,8 @@ ldlPackedSH x = {-# SCC "ldlPacked" #-} LDL m p | |||
376 | (m,p) = ldlPacked' x | 370 | (m,p) = ldlPacked' x |
377 | 371 | ||
378 | -- | Obtains the LDL decomposition of a matrix in a compact data structure suitable for 'ldlSolve'. | 372 | -- | Obtains the LDL decomposition of a matrix in a compact data structure suitable for 'ldlSolve'. |
379 | ldlPacked :: Field t => Her t -> LDL t | 373 | ldlPacked :: Field t => Herm t -> LDL t |
380 | ldlPacked (Her m) = ldlPackedSH m | 374 | ldlPacked (Herm m) = ldlPackedSH m |
381 | 375 | ||
382 | -- | Solution of a linear system (for several right hand sides) from a precomputed LDL factorization obtained by 'ldlPacked'. | 376 | -- | Solution of a linear system (for several right hand sides) from a precomputed LDL factorization obtained by 'ldlPacked'. |
383 | -- | 377 | -- |
@@ -461,26 +455,39 @@ fromList [11.344814282762075,0.17091518882717918,-0.5157294715892575] | |||
461 | 3.000 5.000 6.000 | 455 | 3.000 5.000 6.000 |
462 | 456 | ||
463 | -} | 457 | -} |
464 | eigSH :: Field t => Her t -> (Vector Double, Matrix t) | 458 | eigSH :: Field t => Herm t -> (Vector Double, Matrix t) |
465 | eigSH (Her m) = eigSH' m | 459 | eigSH (Herm m) = eigSH' m |
466 | 460 | ||
467 | -- | Eigenvalues (in descending order) of a complex hermitian or real symmetric matrix. | 461 | -- | Eigenvalues (in descending order) of a complex hermitian or real symmetric matrix. |
468 | eigenvaluesSH :: Field t => Her t -> Vector Double | 462 | eigenvaluesSH :: Field t => Herm t -> Vector Double |
469 | eigenvaluesSH (Her m) = eigenvaluesSH' m | 463 | eigenvaluesSH (Herm m) = eigenvaluesSH' m |
470 | 464 | ||
471 | -------------------------------------------------------------- | 465 | -------------------------------------------------------------- |
472 | 466 | ||
467 | -- | QR decomposition of a matrix in compact form. (The orthogonal matrix is not explicitly formed.) | ||
468 | data QR t = QR (Matrix t) (Vector t) | ||
469 | |||
470 | instance (NFData t, Numeric t) => NFData (QR t) | ||
471 | where | ||
472 | rnf (QR m _) = rnf m | ||
473 | |||
474 | |||
473 | -- | QR factorization. | 475 | -- | QR factorization. |
474 | -- | 476 | -- |
475 | -- If @(q,r) = qr m@ then @m == q \<> r@, where q is unitary and r is upper triangular. | 477 | -- If @(q,r) = qr m@ then @m == q \<> r@, where q is unitary and r is upper triangular. |
476 | qr :: Field t => Matrix t -> (Matrix t, Matrix t) | 478 | qr :: Field t => Matrix t -> (Matrix t, Matrix t) |
477 | qr = {-# SCC "qr" #-} unpackQR . qr' | 479 | qr = {-# SCC "qr" #-} unpackQR . qr' |
478 | 480 | ||
479 | qrRaw m = qr' m | 481 | -- | Compute the QR decomposition of a matrix in compact form. |
482 | qrRaw :: Field t => Matrix t -> QR t | ||
483 | qrRaw m = QR x v | ||
484 | where | ||
485 | (x,v) = qr' m | ||
480 | 486 | ||
481 | {- | generate a matrix with k orthogonal columns from the output of qrRaw | 487 | -- | generate a matrix with k orthogonal columns from the compact QR decomposition obtained by 'qrRaw'. |
482 | -} | 488 | -- |
483 | qrgr n (a,t) | 489 | qrgr :: Field t => Int -> QR t -> Matrix t |
490 | qrgr n (QR a t) | ||
484 | | dim t > min (cols a) (rows a) || n < 0 || n > dim t = error "qrgr expects k <= min(rows,cols)" | 491 | | dim t > min (cols a) (rows a) || n < 0 || n > dim t = error "qrgr expects k <= min(rows,cols)" |
485 | | otherwise = qrgr' n (a,t) | 492 | | otherwise = qrgr' n (a,t) |
486 | 493 | ||
@@ -525,12 +532,12 @@ cholSH = cholSH' | |||
525 | -- | Cholesky factorization of a positive definite hermitian or symmetric matrix. | 532 | -- | Cholesky factorization of a positive definite hermitian or symmetric matrix. |
526 | -- | 533 | -- |
527 | -- If @c = chol m@ then @c@ is upper triangular and @m == tr c \<> c@. | 534 | -- If @c = chol m@ then @c@ is upper triangular and @m == tr c \<> c@. |
528 | chol :: Field t => Her t -> Matrix t | 535 | chol :: Field t => Herm t -> Matrix t |
529 | chol (Her m) = {-# SCC "chol" #-} cholSH' m | 536 | chol (Herm m) = {-# SCC "chol" #-} cholSH' m |
530 | 537 | ||
531 | -- | Similar to 'chol', but instead of an error (e.g., caused by a matrix not positive definite) it returns 'Nothing'. | 538 | -- | Similar to 'chol', but instead of an error (e.g., caused by a matrix not positive definite) it returns 'Nothing'. |
532 | mbChol :: Field t => Her t -> Maybe (Matrix t) | 539 | mbChol :: Field t => Herm t -> Maybe (Matrix t) |
533 | mbChol (Her m) = {-# SCC "mbChol" #-} mbCholSH' m | 540 | mbChol (Herm m) = {-# SCC "mbChol" #-} mbCholSH' m |
534 | 541 | ||
535 | 542 | ||
536 | 543 | ||
@@ -870,6 +877,8 @@ fixPerm' s = res $ mutable f s0 | |||
870 | triang r c h v = (r><c) [el s t | s<-[0..r-1], t<-[0..c-1]] | 877 | triang r c h v = (r><c) [el s t | s<-[0..r-1], t<-[0..c-1]] |
871 | where el p q = if q-p>=h then v else 1 - v | 878 | where el p q = if q-p>=h then v else 1 - v |
872 | 879 | ||
880 | -- | Compute the explicit LU decomposition from the compact one obtained by 'luPacked'. | ||
881 | luFact :: Numeric t => LU t -> (Matrix t, Matrix t, Matrix t, t) | ||
873 | luFact (LU l_u perm) | 882 | luFact (LU l_u perm) |
874 | | r <= c = (l ,u ,p, s) | 883 | | r <= c = (l ,u ,p, s) |
875 | | otherwise = (l',u',p, s) | 884 | | otherwise = (l',u',p, s) |
@@ -967,10 +976,10 @@ relativeError norm a b = r | |||
967 | -- | Generalized symmetric positive definite eigensystem Av = lBv, | 976 | -- | Generalized symmetric positive definite eigensystem Av = lBv, |
968 | -- for A and B symmetric, B positive definite. | 977 | -- for A and B symmetric, B positive definite. |
969 | geigSH :: Field t | 978 | geigSH :: Field t |
970 | => Her t -- ^ A | 979 | => Herm t -- ^ A |
971 | -> Her t -- ^ B | 980 | -> Herm t -- ^ B |
972 | -> (Vector Double, Matrix t) | 981 | -> (Vector Double, Matrix t) |
973 | geigSH (Her a) (Her b) = geigSH' a b | 982 | geigSH (Herm a) (Herm b) = geigSH' a b |
974 | 983 | ||
975 | geigSH' :: Field t | 984 | geigSH' :: Field t |
976 | => Matrix t -- ^ A | 985 | => Matrix t -- ^ A |
@@ -989,29 +998,33 @@ geigSH' a b = (l,v') | |||
989 | 998 | ||
990 | -- | A matrix that, by construction, it is known to be complex Hermitian or real symmetric. | 999 | -- | A matrix that, by construction, it is known to be complex Hermitian or real symmetric. |
991 | -- | 1000 | -- |
992 | -- It can be created using 'sym', 'xTx', or 'trustSym', and the matrix can be extracted using 'her'. | 1001 | -- It can be created using 'sym', 'mTm', or 'trustSym', and the matrix can be extracted using 'unSym'. |
993 | data Her t = Her (Matrix t) deriving Show | 1002 | newtype Herm t = Herm (Matrix t) deriving Show |
1003 | |||
1004 | instance (NFData t, Numeric t) => NFData (Herm t) | ||
1005 | where | ||
1006 | rnf (Herm m) = rnf m | ||
994 | 1007 | ||
995 | -- | Extract the general matrix from a 'Her' structure, forgetting its symmetric or Hermitian property. | 1008 | -- | Extract the general matrix from a 'Herm' structure, forgetting its symmetric or Hermitian property. |
996 | her :: Her t -> Matrix t | 1009 | unSym :: Herm t -> Matrix t |
997 | her (Her x) = x | 1010 | unSym (Herm x) = x |
998 | 1011 | ||
999 | -- | Compute the complex Hermitian or real symmetric part of a square matrix (@(x + tr x)/2@). | 1012 | -- | Compute the complex Hermitian or real symmetric part of a square matrix (@(x + tr x)/2@). |
1000 | sym :: Field t => Matrix t -> Her t | 1013 | sym :: Field t => Matrix t -> Herm t |
1001 | sym x = Her (scale 0.5 (tr x `add` x)) | 1014 | sym x = Herm (scale 0.5 (tr x `add` x)) |
1002 | 1015 | ||
1003 | -- | Compute the contraction @tr x <> x@ of a general matrix. | 1016 | -- | Compute the contraction @tr x <> x@ of a general matrix. |
1004 | xTx :: Numeric t => Matrix t -> Her t | 1017 | mTm :: Numeric t => Matrix t -> Herm t |
1005 | xTx x = Her (tr x `mXm` x) | 1018 | mTm x = Herm (tr x `mXm` x) |
1006 | 1019 | ||
1007 | instance Field t => Linear t Her where | 1020 | instance Field t => Linear t Herm where |
1008 | scale x (Her m) = Her (scale x m) | 1021 | scale x (Herm m) = Herm (scale x m) |
1009 | 1022 | ||
1010 | instance Field t => Additive (Her t) where | 1023 | instance Field t => Additive (Herm t) where |
1011 | add (Her a) (Her b) = Her (a `add` b) | 1024 | add (Herm a) (Herm b) = Herm (a `add` b) |
1012 | 1025 | ||
1013 | -- | At your own risk, declare that a matrix is complex Hermitian or real symmetric | 1026 | -- | At your own risk, declare that a matrix is complex Hermitian or real symmetric |
1014 | -- for usage in 'chol', 'eigSH', etc. Only a triangular part of the matrix will be used. | 1027 | -- for usage in 'chol', 'eigSH', etc. Only a triangular part of the matrix will be used. |
1015 | trustSym :: Matrix t -> Her t | 1028 | trustSym :: Matrix t -> Herm t |
1016 | trustSym x = (Her x) | 1029 | trustSym x = (Herm x) |
1017 | 1030 | ||
diff --git a/packages/base/src/Internal/Modular.hs b/packages/base/src/Internal/Modular.hs index a3421a8..239c742 100644 --- a/packages/base/src/Internal/Modular.hs +++ b/packages/base/src/Internal/Modular.hs | |||
@@ -371,7 +371,7 @@ test = (ok, info) | |||
371 | 371 | ||
372 | checkLU okf t = norm_Inf $ flatten (l <> u <> p - t) | 372 | checkLU okf t = norm_Inf $ flatten (l <> u <> p - t) |
373 | where | 373 | where |
374 | (l,u,p,_ :: Int) = luFact (LU x' p') | 374 | (l,u,p,_) = luFact (LU x' p') |
375 | where | 375 | where |
376 | (x',p') = mutable (luST okf) t | 376 | (x',p') = mutable (luST okf) t |
377 | 377 | ||
diff --git a/packages/base/src/Internal/Util.hs b/packages/base/src/Internal/Util.hs index 36b7855..cf42961 100644 --- a/packages/base/src/Internal/Util.hs +++ b/packages/base/src/Internal/Util.hs | |||
@@ -462,7 +462,7 @@ null1 :: Matrix R -> Vector R | |||
462 | null1 = last . toColumns . snd . rightSV | 462 | null1 = last . toColumns . snd . rightSV |
463 | 463 | ||
464 | -- | solution of overconstrained homogeneous symmetric linear system | 464 | -- | solution of overconstrained homogeneous symmetric linear system |
465 | null1sym :: Her R -> Vector R | 465 | null1sym :: Herm R -> Vector R |
466 | null1sym = last . toColumns . snd . eigSH | 466 | null1sym = last . toColumns . snd . eigSH |
467 | 467 | ||
468 | -------------------------------------------------------------------------------- | 468 | -------------------------------------------------------------------------------- |