sparse-linear-algebra-0.2.0.5: Numerical computation in native Haskell

Safe HaskellNone
LanguageHaskell2010

Numeric.LinearAlgebra.Sparse

Contents

Synopsis

Sparsify : remove almost-0 elements (|x| < eps)

sparsifySV :: SpVector Double -> SpVector Double Source #

Sparsify an SpVector

Matrix condition number

conditionNumberSM :: SpMatrix Double -> Double Source #

uses the R matrix from the QR factorization

Householder transformation

hhMat :: Num a => a -> SpVector a -> SpMatrix a Source #

hhRefl :: SpVector Double -> SpMatrix Double Source #

a vector x uniquely defines an orthogonal plane; the Householder operator reflects any point v with respect to this plane: v' = (I - 2 x >< x) v

Givens rotation matrix

hypot :: Floating a => a -> a -> a Source #

sign :: (Ord a, Num a) => a -> a Source #

givensCoef :: (Ord a, Floating a) => a -> a -> (a, a, a) Source #

Givens coefficients (using stable algorithm shown in Anderson, Edward (4 December 2000). "Discontinuous Plane Rotations and the Symmetric Eigenvalue Problem". LAPACK Working Note)

givens :: SpMatrix Double -> IxRow -> IxCol -> SpMatrix Double Source #

Givens method, row version: choose other row index i' s.t. i' is : * below the diagonal * corresponding element is nonzero

QR.C1 ) To zero out entry A(i, j) we must find row k such that A(k, j) is non-zero but A has zeros in row k for all columns less than j.

firstNonZeroColumn :: IntMap a -> IxRow -> Bool Source #

Is the kth the first nonzero column in the row?

candidateRows :: IntMap (IntMap a) -> IxRow -> IxCol -> Maybe [Key] Source #

Returns a set of rows {k} that satisfy QR.C1

QR decomposition

qr :: SpMatrix Double -> (SpMatrix Double, SpMatrix Double) Source #

Applies Givens rotation iteratively to zero out sub-diagonal elements

gmats :: SpMatrix Double -> [SpMatrix Double] Source #

Givens matrices in order [G1, G2, .. , G_N ]

Eigenvalue algorithms

All eigenvalues (QR algorithm)

eigsQR :: Int -> SpMatrix Double -> SpVector Double Source #

`eigsQR n mm` performs n iterations of the QR algorithm on matrix mm

One eigenvalue and eigenvector (Rayleigh iteration)

eigRayleigh :: Int -> SpMatrix Double -> (SpVector Double, Double) -> (SpVector Double, Double) Source #

`eigsRayleigh n mm` performs n iterations of the Rayleigh algorithm on matrix mm. Cubic-order convergence, but it requires a mildly educated guess on the initial eigenpair

Householder vector

Householder bidiagonalization

SVD

LU factorization

Iterative linear solvers

BCG

bcgStep :: SpMatrix Double -> BCG -> BCG Source #

one step of BCG

data BCG Source #

Constructors

BCG 

Instances

Eq BCG Source # 

Methods

(==) :: BCG -> BCG -> Bool #

(/=) :: BCG -> BCG -> Bool #

Show BCG Source # 

Methods

showsPrec :: Int -> BCG -> ShowS #

show :: BCG -> String #

showList :: [BCG] -> ShowS #

CGS

cgsStep :: SpMatrix Double -> SpVector Double -> CGS -> CGS Source #

one step of CGS

data CGS Source #

Constructors

CGS 

Fields

Instances

Eq CGS Source # 

Methods

(==) :: CGS -> CGS -> Bool #

(/=) :: CGS -> CGS -> Bool #

Show CGS Source # 

Methods

showsPrec :: Int -> CGS -> ShowS #

show :: CGS -> String #

showList :: [CGS] -> ShowS #

cgs :: SpMatrix Double -> SpVector Double -> SpVector Double -> SpVector Double -> CGS Source #

iterate solver until convergence or until max # of iterations is reached

BiCGSTAB

bicgstab :: SpMatrix Double -> SpVector Double -> SpVector Double -> SpVector Double -> BICGSTAB Source #

iterate solver until convergence or until max # of iterations is reached

Moore-Penrose pseudoinverse

pinv :: SpMatrix Double -> SpVector Double -> SpVector Double Source #

Least-squares approximation of a rectangular system of equaitons. Uses \ for the linear solve

Linear solver interface

linSolve :: LinSolveMethod -> SpMatrix Double -> SpVector Double -> SpVector Double Source #

Linear solve with _deterministic_ starting vector (every component at 0.1)

(<\>) :: SpMatrix Double -> SpVector Double -> SpVector Double Source #

\ : linSolve using the BiCGSTAB method as default

Control primitives for bounded iteration with convergence check

modifyUntil :: MonadState s m => (s -> Bool) -> (s -> s) -> m s Source #

transform state until a condition is met

loopUntilAcc :: Int -> ([t] -> Bool) -> (t -> t) -> t -> t Source #

Keep a moving window buffer (length 2) of state x to assess convergence, stop when either a condition on that list is satisfied or when max # of iterations is reached

modifyInspectN :: MonadState s m => Int -> ([s] -> Bool) -> (s -> s) -> m s Source #

Keep a moving window buffer (length 2) of state x to assess convergence, stop when either a condition on that list is satisfied or when max # of iterations is reached (i.e. same thing as loopUntilAcc but this one runs in the State monad)

meanl :: (Foldable t, Fractional a) => t a -> a Source #

norm2l :: (Foldable t, Functor t, Floating a) => t a -> a Source #

diffSqL :: Floating a => [a] -> a Source #

untilConverged :: MonadState a m => (a -> SpVector Double) -> (a -> a) -> m a Source #

iterate until convergence is verified or we run out of a fixed iteration budget

normDiffConverged :: (Foldable t, Functor t) => (a -> SpVector Double) -> t a -> Bool Source #

convergence check (FIXME)

runAppendN :: ([t] -> Bool) -> (t -> t) -> Int -> t -> [t] Source #

run niter iterations and append the state x to a list xs, stop when either the xs satisfies a predicate q or when the counter reaches 0

runAppendN' :: (t -> t) -> Int -> t -> [t] Source #

", NO convergence check

Random matrices and vectors

randMat :: PrimMonad m => Int -> m (SpMatrix Double) Source #

Dense SpMatrix

randVec :: PrimMonad m => Int -> m (SpVector Double) Source #

Dense SpVector

randSpMat :: Int -> Int -> IO (SpMatrix Double) Source #

Sparse SpMatrix

randSpVec :: Int -> Int -> IO (SpVector Double) Source #

Sparse SpVector

Pretty printing

showNonZero :: (Show a, Num a, Eq a) => a -> String Source #

toDenseRow :: Num a => SpMatrix a -> Key -> [a] Source #

toDenseRowClip :: (Show a, Num a) => SpMatrix a -> Key -> Int -> String Source #

printDenseSM :: (Show t, Num t) => SpMatrix t -> IO () Source #

printDenseSV :: (Show t, Num t) => SpVector t -> IO () Source #

Pretty printer typeclass

class PrintDense a where Source #

Minimal complete definition

prd

Methods

prd :: a -> IO () Source #

Instances

(Show a, Num a) => PrintDense (SpVector a) Source # 

Methods

prd :: SpVector a -> IO () Source #

(Show a, Num a) => PrintDense (SpMatrix a) Source # 

Methods

prd :: SpMatrix a -> IO () Source #

Type synonyms for SpMatrix