-- Hoogle documentation, generated by Haddock
-- See Hoogle, http://www.haskell.org/hoogle/
-- | Numerical computation in native Haskell
--
-- Currently the library provides iterative linear solvers, matrix
-- decompositions, eigenvalue computations and related utilities. Please
-- see README.md for details
@package sparse-linear-algebra
@version 0.2.0.7
module Numeric.LinearAlgebra.Class
class Functor f => Additive f
-- | Ring zero element
zero :: (Additive f, Num a) => f a
-- | Ring +
(^+^) :: (Additive f, Num a) => f a -> f a -> f a
-- | negate the values in a functor
negated :: (Num a, Functor f) => f a -> f a
-- | subtract two Additive objects
(^-^) :: (Additive f, Num a) => f a -> f a -> f a
class Additive f => VectorSpace f
-- | multiplication by a scalar
(.*) :: (VectorSpace f, Num a) => a -> f a -> f a
-- | linear interpolation
lerp :: (VectorSpace f, Num a) => a -> f a -> f a -> f a
class VectorSpace f => Hilbert f
-- | inner product
dot :: (Hilbert f, Num a) => f a -> f a -> a
-- | `hilbertDistSq x y = || x - y ||^2`
hilbertDistSq :: (Hilbert f, Num a) => f a -> f a -> a
class Hilbert f => Normed f
norm :: (Normed f, Floating a, Eq a) => a -> f a -> a
-- | Squared 2-norm
normSq :: (Hilbert f, Num a) => f a -> a
-- | L1 norm
norm1 :: (Foldable t, Num a, Functor t) => t a -> a
-- | Euclidean norm
norm2 :: (Hilbert f, Floating a) => f a -> a
-- | Lp norm (p > 0)
normP :: (Foldable t, Functor t, Floating a) => a -> t a -> a
-- | Infinity-norm
normInfty :: (Foldable t, Ord a) => t a -> a
-- | Normalize w.r.t. p-norm (p finite)
normalize :: (Normed f, Floating a, Eq a) => a -> f a -> f a
-- | Lp inner product (p > 0)
dotLp :: (Set t, Foldable t, Floating a) => a -> t a -> t a -> a
-- | Reciprocal
reciprocal :: (Functor f, Fractional b) => f b -> f b
-- | Scale
scale :: (Num b, Functor f) => b -> f b -> f b
class Additive f => FiniteDim f where type FDSize f :: * where {
type family FDSize f :: *;
}
dim :: FiniteDim f => f a -> FDSize f
-- | unary dimension-checking bracket
withDim :: (FiniteDim f, Show e) => f a -> (FDSize f -> f a -> Bool) -> (f a -> c) -> String -> (f a -> e) -> c
-- | binary dimension-checking bracket
withDim2 :: (FiniteDim f, FiniteDim g, Show e) => f a -> g b -> (FDSize f -> FDSize g -> f a -> g b -> Bool) -> (f a -> g b -> c) -> String -> (f a -> g b -> e) -> c
class Additive f => HasData f a where type HDData f a :: * where {
type family HDData f a :: *;
}
dat :: HasData f a => f a -> HDData f a
class (FiniteDim f, HasData f a) => Sparse f a
spy :: (Sparse f a, Fractional b) => f a -> b
class Functor f => Set f
-- | union binary lift : apply function on _union_ of two Sets
liftU2 :: Set f => (a -> a -> a) -> f a -> f a -> f a
-- | intersection binary lift : apply function on _intersection_ of two
-- Sets
liftI2 :: Set f => (a -> b -> c) -> f a -> f b -> f c
class IxContainer (c :: * -> *) a where type Ix c :: * where {
type family Ix c :: *;
}
ixcLookup :: IxContainer c a => Ix c -> c a -> Maybe a
ixcLookupDefault :: IxContainer c a => a -> Ix c -> c a -> a
ixcFilter :: IxContainer c a => (a -> Bool) -> c a -> c a
ixcIfilter :: IxContainer c a => (Ix c -> a -> Bool) -> c a -> c a
ixcInsert :: IxContainer c a => Ix c -> a -> c a -> c a
ixcFromList :: IxContainer c a => [(Ix c, a)] -> c a
ixcToList :: IxContainer c a => c a -> [(Ix c, a)]
module Numeric.LinearAlgebra.Sparse.IntMap
-- | Insert an element
insertIM2 :: Key -> Key -> a -> IntMap (IntMap a) -> IntMap (IntMap a)
-- | Lookup a key
lookupIM2 :: Key -> Key -> IntMap (IntMap a) -> Maybe a
-- | Ppopulate an IM2 from a list of (row index, column index, value)
fromListIM2 :: Foldable t => t (Key, Key, a) -> IntMap (IntMap a) -> IntMap (IntMap a)
-- | Indexed left fold over an IM2, with general accumulator
ifoldlIM2' :: (Key -> Key -> a -> b -> b) -> b -> IntMap (IntMap a) -> b
-- | Indexed left fold over an IM2
ifoldlIM2 :: (Key -> Key -> t -> IntMap a -> IntMap a) -> IntMap (IntMap t) -> IntMap a
-- | Left fold over an IM2, with general accumulator
foldlIM2 :: (a -> b -> b) -> b -> IntMap (IntMap a) -> b
-- | Inner indices become outer ones and vice versa. No loss of information
-- because both inner and outer IntMaps are nubbed.
transposeIM2 :: IntMap (IntMap a) -> IntMap (IntMap a)
-- | Map over outer IM and filter all inner IM's
ifilterIM2 :: (Key -> Key -> a -> Bool) -> IntMap (IntMap a) -> IntMap (IntMap a)
-- | Specialized filtering : keep only sub-diagonal elements
filterSubdiag :: IntMap (IntMap a) -> IntMap (IntMap a)
countSubdiagonalNZ :: IntMap (IntMap a) -> Int
-- | List of (row, col) indices of (nonzero) subdiagonal elements
subdiagIndices :: IntMap (IntMap a) -> [(Key, Key)]
rpairs :: (a, [b]) -> [(a, b)]
-- | Map over IM2
mapIM2 :: (a -> b) -> IntMap (IntMap a) -> IntMap (IntMap b)
-- | Indexed map over IM2
imapIM2 :: (Key -> Key -> a -> b) -> IntMap (IntMap a) -> IntMap (IntMap b)
-- | Mapping keys
mapKeysIM2 :: (Key -> Key) -> (Key -> Key) -> IntMap (IntMap a) -> IntMap (IntMap a)
mapColumnIM2 :: (b -> b) -> IntMap (IntMap b) -> Int -> IntMap (IntMap b)
instance Numeric.LinearAlgebra.Class.Set Data.IntMap.Base.IntMap
instance Numeric.LinearAlgebra.Class.Additive Data.IntMap.Base.IntMap
instance Numeric.LinearAlgebra.Class.VectorSpace Data.IntMap.Base.IntMap
instance Numeric.LinearAlgebra.Class.Hilbert Data.IntMap.Base.IntMap
instance Numeric.LinearAlgebra.Class.Normed Data.IntMap.Base.IntMap
module Numeric.Eps
-- | eps = 1e-8
eps :: Double
-- | Rounding rule
almostZero :: Double -> Bool
-- | Rounding rule
almostOne :: Double -> Bool
isNz :: Double -> Bool
withDefault :: (t -> Bool) -> t -> t -> t
roundZero :: Double -> Double
roundOne :: Double -> Double
with2Defaults :: (t -> Bool) -> (t -> Bool) -> t -> t -> t -> t
-- | Round to respectively 0 or 1 within some predefined numerical
-- precision eps
roundZeroOne :: Double -> Double
module Data.Sparse.Utils
-- | Componentwise tuple operations TODO : use semilattice properties
-- instead
maxTup :: Ord t => (t, t) -> (t, t) -> (t, t)
-- | Componentwise tuple operations TODO : use semilattice properties
-- instead
minTup :: Ord t => (t, t) -> (t, t) -> (t, t)
-- | integer-indexed ziplist
denseIxArray :: [b] -> [(Int, b)]
-- | ", 2d arrays
denseIxArray2 :: Int -> [c] -> [(Int, Int, c)]
-- | foldr over the results of a fmap
foldrMap :: (Foldable t, Functor t) => (a -> b) -> (b -> c -> c) -> c -> t a -> c
-- | strict left fold
foldlStrict :: (a -> b -> a) -> a -> [b] -> a
-- | indexed right fold
ifoldr :: Num i => (a -> b -> b) -> b -> (i -> c -> d -> a) -> c -> [d] -> b
type LB = Int
type UB = Int
inBounds :: LB -> UB -> Int -> Bool
inBounds2 :: (LB, UB) -> (Int, Int) -> Bool
inBounds0 :: UB -> Int -> Bool
inBounds02 :: (UB, UB) -> (Int, Int) -> Bool
module Data.Sparse.Types
type Rows = Int
type Cols = Int
type IxRow = Int
type IxCol = Int
module Data.Sparse.SpVector
data SpVector a
SV :: Int -> IntMap a -> SpVector a
[svDim] :: SpVector a -> Int
[svData] :: SpVector a -> IntMap a
-- | SpVector sparsity
spySV :: Fractional b => SpVector a -> b
-- | Number of nonzeros
nzSV :: SpVector a -> Int
sizeStrSV :: SpVector a -> String
-- | empty sparse vector (length n, no entries)
zeroSV :: Int -> SpVector a
-- | singleton sparse vector (length 1)
singletonSV :: a -> SpVector a
-- | create a sparse vector from an association list while discarding all
-- zero entries
mkSpVector :: (Num a, Eq a) => Int -> IntMap a -> SpVector a
-- | ", from logically dense array (consecutive indices)
mkSpVectorD :: (Num a, Eq a) => Int -> [a] -> SpVector a
mkSpVector1 :: Int -> IntMap a -> SpVector a
-- | Create new sparse vector, assumin 0-based, contiguous indexing
fromListDenseSV :: Int -> [a] -> SpVector a
-- | one-hot encoding : `oneHotSV n k` produces a SpVector of length n
-- having 1 at the k-th position
oneHotSVU :: Num a => Int -> IxRow -> SpVector a
oneHotSV :: Num a => Int -> IxRow -> SpVector a
-- | DENSE vector of `1`s
onesSV :: Num a => Int -> SpVector a
-- | DENSE vector of `0`s
zerosSV :: Num a => Int -> SpVector a
-- | insert element x at index i in a preexisting
-- SpVector
insertSpVector :: Int -> a -> SpVector a -> SpVector a
fromListSV :: Int -> [(Int, a)] -> SpVector a
toListSV :: SpVector a -> [(Key, a)]
-- | To dense list (default = 0)
toDenseListSV :: Num b => SpVector b -> [b]
-- | Lookup an index in a SpVector
lookupSV :: Key -> SpVector a -> Maybe a
-- | Lookup an index, return a default value if lookup fails
lookupDefaultSV :: a -> Key -> SpVector a -> a
-- | Lookup an index in a SpVector, returns 0 if lookup fails
lookupDenseSV :: Num a => Key -> SpVector a -> a
-- | Tail elements
tailSV :: SpVector a -> SpVector a
-- | Head element
headSV :: Num a => SpVector a -> a
-- | Concatenate two sparse vectors
concatSV :: SpVector a -> SpVector a -> SpVector a
-- | Filter
filterSV :: (a -> Bool) -> SpVector a -> SpVector a
-- | Indexed filter
ifilterSV :: (Int -> a -> Bool) -> SpVector a -> SpVector a
-- | Generate an arbitrary (not random) vector u such that `v dot
-- u = 0`
orthogonalSV :: Fractional a => SpVector a -> SpVector a
instance GHC.Classes.Eq a => GHC.Classes.Eq (Data.Sparse.SpVector.SpVector a)
instance GHC.Base.Functor Data.Sparse.SpVector.SpVector
instance Numeric.LinearAlgebra.Class.Set Data.Sparse.SpVector.SpVector
instance Data.Foldable.Foldable Data.Sparse.SpVector.SpVector
instance Numeric.LinearAlgebra.Class.Additive Data.Sparse.SpVector.SpVector
instance Numeric.LinearAlgebra.Class.VectorSpace Data.Sparse.SpVector.SpVector
instance Numeric.LinearAlgebra.Class.FiniteDim Data.Sparse.SpVector.SpVector
instance Numeric.LinearAlgebra.Class.HasData Data.Sparse.SpVector.SpVector a
instance Numeric.LinearAlgebra.Class.Sparse Data.Sparse.SpVector.SpVector a
instance Numeric.LinearAlgebra.Class.Hilbert Data.Sparse.SpVector.SpVector
instance Numeric.LinearAlgebra.Class.Normed Data.Sparse.SpVector.SpVector
instance GHC.Show.Show a => GHC.Show.Show (Data.Sparse.SpVector.SpVector a)
module Data.Sparse.SpMatrix
data SpMatrix a
SM :: (Rows, Cols) -> IntMap (IntMap a) -> SpMatrix a
[smDim] :: SpMatrix a -> (Rows, Cols)
[smData] :: SpMatrix a -> IntMap (IntMap a)
sizeStr :: SpMatrix a -> String
-- | `zeroSM m n` : Empty SpMatrix of size (m, n)
zeroSM :: Rows -> Cols -> SpMatrix a
mkDiagonal :: Int -> [a] -> SpMatrix a
-- | `eye n` : identity matrix of rank n
eye :: Num a => Int -> SpMatrix a
-- | Permutation matrix from a (possibly incomplete) list of row swaps
-- starting from row 0 e.g. `permutationSM 5 [1,3]` first swaps rows (0,
-- 1) and then rows (1, 3) :
--
--
-- - 0,1,0,0,0
-- - 0,0,0,1,0
-- - 0,0,1,0,0
-- - 1,0,0,0,0
-- - 0,0,0,0,1
--
permutationSM :: Num a => Int -> [IxRow] -> SpMatrix a
-- | Permutation matrix from a (possibly incomplete) list of row pair swaps
-- e.g. `permutPairs 5 [(2,4)]` swaps rows (2, 4) :
--
--
-- - 1,0,0,0,0
-- - 0,1,0,0,0
-- - 0,0,0,0,1
-- - 0,0,0,1,0
-- - 0,0,1,0,0
--
permutPairsSM :: Num a => Int -> [(IxRow, IxRow)] -> SpMatrix a
-- | `mkSubDiagonal n o xx` creates a square SpMatrix of size n
-- with xx on the oth subdiagonal
mkSubDiagonal :: Int -> Int -> [a] -> SpMatrix a
-- | Insert an element in a preexisting Spmatrix at the specified indices
insertSpMatrix :: IxRow -> IxCol -> a -> SpMatrix a -> SpMatrix a
-- | Add to existing SpMatrix using data from list (row, col, value)
fromListSM' :: Foldable t => t (IxRow, IxCol, a) -> SpMatrix a -> SpMatrix a
-- | Create new SpMatrix using data from list (row, col, value)
fromListSM :: Foldable t => (Int, Int) -> t (IxRow, IxCol, a) -> SpMatrix a
-- | Create new SpMatrix assuming contiguous, 0-based indexing of elements
fromListDenseSM :: Int -> [a] -> SpMatrix a
-- | Populate list with SpMatrix contents and populate missing entries with
-- 0
toDenseListSM :: Num t => SpMatrix t -> [(IxRow, IxCol, t)]
lookupSM :: SpMatrix a -> IxRow -> IxCol -> Maybe a
-- | Looks up an element in the matrix with a default (if the element is
-- not found, zero is returned)
lookupWD_SM :: Num a => SpMatrix a -> (IxRow, IxCol) -> a
-- | Zero-default lookup, infix form (no bound checking)
--
-- Looks up an element in the matrix with a default (if the element is
-- not found, zero is returned)
(@@!) :: Num a => SpMatrix a -> (IxRow, IxCol) -> a
-- | Zero-default lookup, infix form ("safe" : throws exception if lookup
-- is outside matrix bounds)
--
-- Looks up an element in the matrix with a default (if the element is
-- not found, zero is returned)
(@@) :: Num a => SpMatrix a -> (IxRow, IxCol) -> a
lookupWD_IM :: Num a => IntMap (IntMap a) -> (IxRow, IxCol) -> a
-- | Extract a submatrix given the specified index bounds, rebalancing keys
-- with the two supplied functions
extractSubmatrixSM :: (Key -> Key) -> (Key -> Key) -> SpMatrix a -> (IxRow, IxRow) -> (IxCol, IxCol) -> SpMatrix a
-- | Extract a submatrix given the specified index bounds NB : subtracts
-- (i1, j1) from the indices
extractSubmatrixRebalanceKeys :: SpMatrix a -> (IxRow, IxRow) -> (IxCol, IxCol) -> SpMatrix a
-- | Extract a submatrix given the specified index bounds NB : submatrix
-- indices are _preserved_
extractSubmatrix :: SpMatrix a -> (IxRow, IxRow) -> (IxCol, IxCol) -> SpMatrix a
extractRowSM :: SpMatrix a -> IxRow -> SpMatrix a
-- | Extract column within a row range
extractSubRowSM :: SpMatrix a -> IxRow -> (IxCol, IxCol) -> SpMatrix a
-- | Extract column within a row range, rebalance keys
extractSubRowSM_RK :: SpMatrix a -> IxRow -> (IxCol, IxCol) -> SpMatrix a
-- | Extract all column
extractColSM :: SpMatrix a -> IxCol -> SpMatrix a
-- | Extract column within a row range
extractSubColSM :: SpMatrix a -> IxCol -> (IxRow, IxRow) -> SpMatrix a
-- | Extract column within a row range, rebalance keys
extractSubColSM_RK :: SpMatrix a -> IxCol -> (IxRow, IxRow) -> SpMatrix a
-- | Are the supplied indices within matrix bounds?
isValidIxSM :: SpMatrix a -> (Int, Int) -> Bool
-- | Is the matrix square?
isSquareSM :: SpMatrix a -> Bool
-- | Is the matrix diagonal?
isDiagonalSM :: SpMatrix a -> Bool
-- | is the matrix orthogonal? i.e. Q^t ## Q == I
isOrthogonalSM :: SpMatrix Double -> Bool
-- | Data in internal representation (do not export)
immSM :: SpMatrix t -> IntMap (IntMap t)
-- | (Number of rows, Number of columns)
dimSM :: SpMatrix t -> (Rows, Cols)
-- | Number of rows times number of columns
nelSM :: SpMatrix t -> Int
-- | Number of rows
nrows :: SpMatrix a -> Rows
-- | Number of columns
ncols :: SpMatrix a -> Cols
data SMInfo
SMInfo :: Int -> Double -> SMInfo
[smNz] :: SMInfo -> Int
[smSpy] :: SMInfo -> Double
infoSM :: SpMatrix a -> SMInfo
nzSM :: SpMatrix a -> Int
spySM :: Fractional b => SpMatrix a -> b
nzRow :: SpMatrix a -> Key -> Int
bwMinSM :: SpMatrix a -> Int
bwMaxSM :: SpMatrix a -> Int
bwBoundsSM :: SpMatrix a -> (Int, Int)
-- | Vertical stacking
vertStackSM :: SpMatrix a -> SpMatrix a -> SpMatrix a
-- | Vertical stacking
(-=-) :: SpMatrix a -> SpMatrix a -> SpMatrix a
-- | Horizontal stacking
horizStackSM :: SpMatrix a -> SpMatrix a -> SpMatrix a
-- | Horizontal stacking
(-||-) :: SpMatrix a -> SpMatrix a -> SpMatrix a
-- | Left fold over SpMatrix
foldlSM :: (a -> b -> b) -> b -> SpMatrix a -> b
-- | Indexed left fold over SpMatrix
ifoldlSM :: (Key -> Key -> a -> b -> b) -> b -> SpMatrix a -> b
-- | Count sub-diagonal nonzeros
countSubdiagonalNZSM :: SpMatrix a -> Int
-- | Filter the index subset that lies below the diagonal (used in the QR
-- decomposition, for example)
subdiagIndicesSM :: SpMatrix a -> [(IxRow, IxCol)]
sparsifyIM2 :: IntMap (IntMap Double) -> IntMap (IntMap Double)
-- | Sparsify an SpMatrix
sparsifySM :: SpMatrix Double -> SpMatrix Double
-- | Round almost-0 and almost-1 to 0 and 1 respectively
roundZeroOneSM :: SpMatrix Double -> SpMatrix Double
-- | swap two rows of a SpMatrix (bounds not checked)
swapRows :: IxRow -> IxRow -> SpMatrix a -> SpMatrix a
-- | swap two rows of a SpMatrix (bounds checked)
swapRowsSafe :: IxRow -> IxRow -> SpMatrix a -> SpMatrix a
-- | transposeSM, (#^) : Matrix transpose
transposeSM :: SpMatrix a -> SpMatrix a
-- | transposeSM, (#^) : Matrix transpose
(#^) :: SpMatrix a -> SpMatrix a
matScale :: Num a => a -> SpMatrix a -> SpMatrix a
normFrobenius :: SpMatrix Double -> Double
matMat :: Num a => SpMatrix a -> SpMatrix a -> SpMatrix a
(##) :: Num a => SpMatrix a -> SpMatrix a -> SpMatrix a
-- | Removes all elements x for which `| x | <= eps`)
matMatSparsified :: SpMatrix Double -> SpMatrix Double -> SpMatrix Double
-- | Removes all elements x for which `| x | <= eps`)
(#~#) :: SpMatrix Double -> SpMatrix Double -> SpMatrix Double
-- | A^T B
(#^#) :: SpMatrix Double -> SpMatrix Double -> SpMatrix Double
-- | A B^T
(##^) :: SpMatrix Double -> SpMatrix Double -> SpMatrix Double
-- | Contract two matrices A and B up to an index n, i.e. summing
-- over repeated indices: Aij Bjk , for j in [0 .. n]
contractSub :: Num a => SpMatrix a -> SpMatrix a -> IxRow -> IxCol -> Int -> a
instance GHC.Show.Show Data.Sparse.SpMatrix.SMInfo
instance GHC.Classes.Eq Data.Sparse.SpMatrix.SMInfo
instance GHC.Classes.Eq a => GHC.Classes.Eq (Data.Sparse.SpMatrix.SpMatrix a)
instance GHC.Show.Show a => GHC.Show.Show (Data.Sparse.SpMatrix.SpMatrix a)
instance GHC.Base.Functor Data.Sparse.SpMatrix.SpMatrix
instance Numeric.LinearAlgebra.Class.Set Data.Sparse.SpMatrix.SpMatrix
instance Numeric.LinearAlgebra.Class.Additive Data.Sparse.SpMatrix.SpMatrix
instance Numeric.LinearAlgebra.Class.FiniteDim Data.Sparse.SpMatrix.SpMatrix
instance Numeric.LinearAlgebra.Class.HasData Data.Sparse.SpMatrix.SpMatrix a
instance Numeric.LinearAlgebra.Class.Sparse Data.Sparse.SpMatrix.SpMatrix a
module Data.Sparse.Common
-- | Insert row , using the provided row index transformation function
insertRowWith :: (IxCol -> IxCol) -> SpMatrix a -> SpVector a -> Key -> SpMatrix a
-- | Insert row
insertRow :: SpMatrix a -> SpVector a -> Key -> SpMatrix a
-- | Insert column, using the provided row index transformation function
insertColWith :: (IxRow -> IxRow) -> SpMatrix a -> SpVector a -> IxCol -> SpMatrix a
-- | Insert column
insertCol :: SpMatrix a -> SpVector a -> IxCol -> SpMatrix a
outerProdSV :: Num a => SpVector a -> SpVector a -> SpMatrix a
(><) :: Num a => SpVector a -> SpVector a -> SpMatrix a
-- | Demote (n x 1) or (1 x n) SpMatrix to SpVector
toSV :: SpMatrix a -> SpVector a
-- | promote a SV to SM
svToSM :: SpVector a -> SpMatrix a
-- | Extract jth column
extractCol :: SpMatrix a -> IxCol -> SpVector a
-- | Extract ith row
extractRow :: SpMatrix a -> IxRow -> SpVector a
-- | Generic extraction function
extractVectorDenseWith :: Num a => (Int -> (IxRow, IxCol)) -> SpMatrix a -> SpVector a
-- | Extract ith row (dense)
extractRowDense :: Num a => SpMatrix a -> IxRow -> SpVector a
-- | Extract jth column
extractColDense :: Num a => SpMatrix a -> IxCol -> SpVector a
-- | Extract the diagonal
extractDiagDense :: Num a => SpMatrix a -> SpVector a
-- | extract row interval
extractSubRow :: SpMatrix a -> IxRow -> (IxCol, IxCol) -> SpVector a
-- | extract column interval
extractSubCol :: SpMatrix a -> IxCol -> (IxRow, IxRow) -> SpVector a
-- | extract row interval, rebalance keys by subtracting lowest one
extractSubRow_RK :: SpMatrix a -> IxRow -> (IxCol, IxCol) -> SpVector a
-- | extract column interval, rebalance keys by subtracting lowest one
extractSubCol_RK :: SpMatrix a -> IxCol -> (IxRow, IxRow) -> SpVector a
-- | Matrix-on-vector
matVec :: Num a => SpMatrix a -> SpVector a -> SpVector a
-- | Matrix-on-vector
(#>) :: Num a => SpMatrix a -> SpVector a -> SpVector a
-- | Vector-on-matrix (FIXME : transposes matrix: more costly than
-- matVec, I think)
vecMat :: Num a => SpVector a -> SpMatrix a -> SpVector a
-- | Vector-on-matrix (FIXME : transposes matrix: more costly than
-- matVec, I think)
(<#) :: Num a => SpVector a -> SpMatrix a -> SpVector a
prd :: PrintDense a => a -> IO ()
instance (GHC.Show.Show a, GHC.Num.Num a) => Data.Sparse.Common.PrintDense (Data.Sparse.SpVector.SpVector a)
instance (GHC.Show.Show a, GHC.Num.Num a) => Data.Sparse.Common.PrintDense (Data.Sparse.SpMatrix.SpMatrix a)
module Numeric.LinearAlgebra.Sparse
-- | Sparsify an SpVector
sparsifySV :: SpVector Double -> SpVector Double
-- | uses the R matrix from the QR factorization
conditionNumberSM :: SpMatrix Double -> Double
hhMat :: Num a => a -> SpVector a -> SpMatrix a
-- | 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
hhRefl :: SpVector Double -> SpMatrix Double
-- | 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.
givens :: SpMatrix Double -> IxRow -> IxCol -> SpMatrix Double
-- | `eigsQR n mm` performs n iterations of the QR algorithm on
-- matrix mm, and returns a SpVector containing all eigenvalues
eigsQR :: Int -> SpMatrix Double -> SpVector Double
-- | `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
eigRayleigh :: Int -> SpMatrix Double -> (SpVector Double, Double) -> (SpVector Double, Double)
cgne :: SpMatrix Double -> SpVector Double -> SpVector Double -> CGNE
tfqmr :: SpMatrix Double -> SpVector Double -> SpVector Double -> TFQMR
-- | iterate solver until convergence or until max # of iterations is
-- reached
bicgstab :: SpMatrix Double -> SpVector Double -> SpVector Double -> SpVector Double -> BICGSTAB
-- | iterate solver until convergence or until max # of iterations is
-- reached
cgs :: SpMatrix Double -> SpVector Double -> SpVector Double -> SpVector Double -> CGS
bcg :: SpMatrix Double -> SpVector Double -> SpVector Double -> BCG
_xCgne :: CGNE -> SpVector Double
_xTfq :: TFQMR -> SpVector Double
_xBicgstab :: BICGSTAB -> SpVector Double
_x :: CGS -> SpVector Double
_xBcg :: BCG -> SpVector Double
-- | one step of CGS
cgsStep :: SpMatrix Double -> SpVector Double -> CGS -> CGS
-- | one step of BiCGSTAB
bicgstabStep :: SpMatrix Double -> SpVector Double -> BICGSTAB -> BICGSTAB
data CGNE
data TFQMR
data BICGSTAB
data CGS
data BCG
-- | Linear solve with _deterministic_ starting vector (every component at
-- 0.1)
linSolve :: LinSolveMethod -> SpMatrix Double -> SpVector Double -> SpVector Double
data LinSolveMethod
-- | \ : linSolve using the BiCGSTAB method as default
(<\>) :: SpMatrix Double -> SpVector Double -> SpVector Double
-- | Dense SpMatrix
randMat :: PrimMonad m => Int -> m (SpMatrix Double)
-- | Dense SpVector
randVec :: PrimMonad m => Int -> m (SpVector Double)
-- | Sparse SpMatrix
randSpMat :: Int -> Int -> IO (SpMatrix Double)
-- | Sparse SpVector
randSpVec :: Int -> Int -> IO (SpVector Double)
-- | 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)
modifyInspectN :: MonadState s m => Int -> ([s] -> Bool) -> (s -> s) -> m s
-- | ", NO convergence check
runAppendN' :: (t -> t) -> Int -> t -> [t]
diffSqL :: Floating a => [a] -> a
-- | Applies Givens rotation iteratively to zero out sub-diagonal elements
qr :: SpMatrix Double -> (SpMatrix Double, SpMatrix Double)
-- | LU factors
lu :: SpMatrix Double -> (SpMatrix Double, SpMatrix Double)
instance GHC.Show.Show Numeric.LinearAlgebra.Sparse.LinSolveMethod
instance GHC.Classes.Eq Numeric.LinearAlgebra.Sparse.LinSolveMethod
instance GHC.Classes.Eq Numeric.LinearAlgebra.Sparse.BICGSTAB
instance GHC.Classes.Eq Numeric.LinearAlgebra.Sparse.CGS
instance GHC.Classes.Eq Numeric.LinearAlgebra.Sparse.BCG
instance GHC.Classes.Eq Numeric.LinearAlgebra.Sparse.TFQMR
instance GHC.Classes.Eq Numeric.LinearAlgebra.Sparse.CGNE
instance GHC.Show.Show Numeric.LinearAlgebra.Sparse.BCG
instance GHC.Show.Show Numeric.LinearAlgebra.Sparse.CGS
instance GHC.Show.Show Numeric.LinearAlgebra.Sparse.BICGSTAB