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

Safe HaskellNone
LanguageHaskell2010

Numeric.LinearAlgebra.Sparse

Contents

Synopsis

Matrix factorizations

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

Applies Givens rotation iteratively to zero out sub-diagonal elements

Incomplete LU

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

used for Incomplete LU : remove entries in m corresponding to zero entries in m2

Condition number

conditionNumberSM :: SpMatrix Double -> Double Source #

uses the R matrix from the QR factorization

Householder reflection

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

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.

Eigensolvers

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

`eigsQR n mm` performs n iterations of the QR algorithm on matrix mm, and returns a SpVector containing all eigenvalues

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

`eigsRayleigh n mm` performs n iterations of the Rayleigh algorithm on matrix mm and returns the eigenpair closest to the initialization. It displays cubic-order convergence, but it also requires an educated guess on the initial eigenpair

Linear solvers

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

Methods

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

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

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

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

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

one step of CGS

data CGNE Source #

Instances

Eq CGNE Source # 

Methods

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

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

data TFQMR Source #

Instances

Eq TFQMR Source # 

Methods

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

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

data CGS Source #

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 #

data BCG Source #

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 #

Matrix partitioning

diagPartitions :: SpMatrix a -> (SpMatrix a, SpMatrix a, SpMatrix a) Source #

Partition a matrix into strictly subdiagonal, diagonal and strictly superdiagonal parts

Random arrays

Random matrices and vectors

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

Dense SpMatrix

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

Dense SpVector

Sparse "

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

Sparse SpMatrix

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

Sparse SpVector

Sparsify data

sparsifySV :: SpVector Double -> SpVector Double Source #

Sparsify an SpVector

Iteration combinators

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)

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

", NO convergence check

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