-- 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.1.1 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 one :: (Additive f, Num a) => f a (^*^) :: (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 Sparse c a => SpContainer c a where type ScIx c :: * where { type family ScIx c :: *; } scInsert :: SpContainer c a => ScIx c -> a -> c a -> c a scLookup :: SpContainer c a => c a -> ScIx c -> Maybe a (@@) :: SpContainer c a => c a -> ScIx c -> a -- | Testing for values "near" zero module Numeric.Eps -- | Provides a test to see if a quantity is near zero. -- --
--   >>> nearZero (1e-11 :: Double)
--   False
--   
-- --
--   >>> nearZero (1e-17 :: Double)
--   True
--   
-- --
--   >>> nearZero (1e-5 :: Float)
--   False
--   
-- --
--   >>> nearZero (1e-7 :: Float)
--   True
--   
class Num a => Epsilon a -- | Determine if a quantity is near zero. nearZero :: Epsilon a => a -> Bool -- | Determine if a quantity is near zero. nearZero :: Epsilon a => a -> Bool -- | Rounding rule isNz :: Epsilon a => a -> Bool roundZero :: Epsilon a => a -> a roundOne :: Epsilon a => a -> a -- | Round to respectively 0 or 1 roundZeroOne :: Epsilon a => a -> a instance Numeric.Eps.Epsilon GHC.Types.Float instance Numeric.Eps.Epsilon GHC.Types.Double instance Numeric.Eps.Epsilon Foreign.C.Types.CFloat instance Numeric.Eps.Epsilon Foreign.C.Types.CDouble module Data.Sparse.Utils -- | Wrap a function with a null check, returning in Maybe harness :: (t -> Bool) -> (t -> a) -> t -> Maybe a -- | 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 head' :: Vector a -> Maybe a tail' :: Vector a -> Maybe (Vector a) module Data.Sparse.Types type Rows = Int type Cols = Int type IxRow = Int type IxCol = Int module Data.Sparse.IntMap2.IntMap2 -- | Insert an element insertIM2 :: Key -> Key -> a -> IntMap (IntMap a) -> IntMap (IntMap a) -- | Lookup a key lookupIM2 :: Key -> Key -> IntMap (IntMap a) -> Maybe a -- | Lookup with default 0 lookupWD_IM :: Num a => IntMap (IntMap a) -> (IxRow, IxCol) -> a -- | Populate 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 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) : -- -- 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 and 4 : -- -- 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 -- | Indexed filtering function filterSM :: (Key -> Key -> a -> Bool) -> SpMatrix a -> SpMatrix a -- | Diagonal, subdiagonal, superdiagonal partitions of a SpMatrix (useful -- for writing preconditioners) extractDiag :: SpMatrix a -> SpMatrix a -- | Diagonal, subdiagonal, superdiagonal partitions of a SpMatrix (useful -- for writing preconditioners) extractSuperDiag :: SpMatrix a -> SpMatrix a -- | Diagonal, subdiagonal, superdiagonal partitions of a SpMatrix (useful -- for writing preconditioners) extractSubDiag :: SpMatrix a -> SpMatrix 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 -- | Extract whole 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 lower/upper triangular? isLowerTriSM :: Eq a => SpMatrix a -> Bool -- | Is the matrix lower/upper triangular? isUpperTriSM :: Eq a => SpMatrix a -> Bool -- | Is the matrix orthogonal? i.e. Q^t ## Q == I isOrthogonalSM :: (Eq a, Epsilon a) => SpMatrix a -> 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 ifilterSM :: (Key -> Key -> a -> Bool) -> 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 :: Epsilon a => IntMap (IntMap a) -> IntMap (IntMap a) -- | Sparsify an SpMatrix sparsifySM :: Epsilon a => SpMatrix a -> SpMatrix a -- | Round almost-0 and almost-1 to 0 and 1 respectively roundZeroOneSM :: Epsilon a => SpMatrix a -> SpMatrix a -- | 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 matScale :: Num a => a -> SpMatrix a -> SpMatrix a normFrobenius :: Floating a => SpMatrix a -> a 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 :: Epsilon a => SpMatrix a -> SpMatrix a -> SpMatrix a -- | Removes all elements x for which `| x | <= eps`) (#~#) :: Epsilon a => SpMatrix a -> SpMatrix a -> SpMatrix a -- | A^T B (#^#) :: Epsilon a => SpMatrix a -> SpMatrix a -> SpMatrix a -- | A B^T (##^) :: Epsilon a => SpMatrix a -> SpMatrix a -> SpMatrix a -- | Contract row i of A with column j of 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 instance GHC.Num.Num a => Numeric.LinearAlgebra.Class.SpContainer Data.Sparse.SpMatrix.SpMatrix a 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 :: Epsilon a => Int -> IntMap a -> SpVector a -- | ", from logically dense array (consecutive indices) mkSpVectorD :: Epsilon 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 -- | Map a function over a range of indices and filter the result (indices -- and values) to fit in a n-long SpVector spVectorDenseIx :: Epsilon a => (Int -> a) -> UB -> [Int] -> SpVector a -- | ", using just the integer bounds of the interval spVectorDenseLoHi :: Epsilon a => (Int -> a) -> UB -> Int -> Int -> 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 -- | Populate a SpVector with the contents of a Vector. fromVector :: Vector a -> SpVector a -- | Populate a Vector with the entries of a SpVector, discarding the -- indices (NB: loses sparsity information). toVector :: SpVector a -> Vector a -- | toVectorDense :: Num a => SpVector a -> Vector 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] -- | Indexed fold over SpVector ifoldSV :: (Key -> a -> b -> b) -> b -> SpVector a -> 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 -- | Keep the first n components of the SpVector (like take for -- lists) takeSV :: Int -> SpVector a -> SpVector a -- | Discard the first n components of the SpVector and rebalance the keys -- (like drop for lists) -- -- Keep the first n components of the SpVector (like take for -- lists) dropSV :: Int -> SpVector a -> SpVector 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 GHC.Num.Num a => Numeric.LinearAlgebra.Class.SpContainer 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.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 -- | Fill the diagonal of a SpMatrix with the components of a SpVector diagonalSM :: SpVector a -> SpMatrix a -- | Outer product (all-with-all matrix) outerProdSV :: Num a => SpVector a -> SpVector a -> SpMatrix a -- | Outer product (all-with-all matrix) (><) :: 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 -- | Lookup a row in a SpMatrix; returns an SpVector with the row, if this -- is non-empty lookupRowSM :: SpMatrix a -> IxRow -> Maybe (SpVector 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 (all entries between columns j1 and j2, INCLUDED, -- are returned) extractSubRow :: SpMatrix a -> IxRow -> (IxCol, -- IxCol) -> SpVector a extractSubRow m i (j1, j2) = case lookupRowSM -- m i of Nothing -> zeroSV (ncols m) Just rv -> ifilterSV (j _ -- -> j >= j1 && j <= j2) rv -- -- ", returning in Maybe extractSubRow :: SpMatrix a -> IxRow -> -- (Int, Int) -> Maybe (SpVector a) extractSubRow m i (j1, j2) = -- resizeSV (j2 - j1) . ifilterSV (j _ -> j >= j1 && j -- j2) <$ lookupRowSM m i -- -- Extract an interval of SpVector components, changing accordingly the -- resulting SpVector size. Keys are _not_ rebalanced, i.e. components -- are still labeled according with respect to the source matrix. extractSubRow :: SpMatrix a -> IxRow -> (Int, Int) -> 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 -- | Pack a V.Vector of SpVectors as columns of an SpMatrix fromCols :: Vector (SpVector a) -> SpMatrix 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 -- | Given a matrix A, returns a pair of matrices (Q, R) such that Q R = A, -- Q is orthogonal and R is upper triangular. Applies Givens rotation -- iteratively to zero out sub-diagonal elements qr :: (Epsilon a, Floating a, Real a) => SpMatrix a -> (SpMatrix a, SpMatrix a) -- | Given a matrix A, returns a pair of matrices (L, U) such that L U = A lu :: (Epsilon a, Fractional a, Real a) => SpMatrix a -> (SpMatrix a, SpMatrix a) -- | Given a positive semidefinite matrix A, returns a lower-triangular -- matrix L such that L L^T = A chol :: (Epsilon a, Real a, Floating a) => SpMatrix a -> SpMatrix a -- | uses the R matrix from the QR factorization conditionNumberSM :: (Epsilon a, RealFloat a) => SpMatrix a -> a 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 :: Num a => SpVector a -> SpMatrix a -- | 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. -- -- NB: the current version is quite inefficient in that: 1. the Givens' -- matrix G_i is different from Identity only in 4 entries 2. at -- each iteration i we multiply G_i by the previous -- partial result M. Since this corresponds to a rotation, and -- the givensCoef function already computes the value of the -- resulting non-zero component (output r), `G_i ## M` can be -- simplified by just changing two entries of M (i.e. zeroing -- one out and changing the other into r). givens :: (Floating a, Epsilon a, Ord a) => SpMatrix a -> IxRow -> IxCol -> SpMatrix a -- | Given a matrix A and a positive integer n, this procedure -- finds the basis of an order n Krylov subspace (as the columns -- of matrix Q), along with an upper Hessenberg matrix H, such that A = -- Q^T H Q. At the i`th iteration, it finds (i + 1) coefficients (the -- i`th column of the Hessenberg matrix H) and the (i + 1)`th Krylov -- vector. arnoldi :: (Floating a, Eq a) => SpMatrix a -> Int -> (SpMatrix a, SpMatrix a) -- | `eigsQR n mm` performs n iterations of the QR algorithm on -- matrix mm, and returns a SpVector containing all eigenvalues eigsQR :: (Epsilon a, Real a, Floating a) => Int -> SpMatrix a -> SpVector a -- | `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 -- -- max # iterations -> SpMatrix Double -- matrix -> (SpVector -- Double, Double) -- initial guess of (eigenvector, eigenvalue) -> -- (SpVector Double, Double) -- final estimate of (eigenvector, -- eigenvalue) eigRayleigh :: Int -> SpMatrix Double -> (SpVector Double, Double) -> (SpVector Double, Double) -- | 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 -- | Direct solver based on a triangular factorization of the system -- matrix. luSolve :: (Fractional a, Eq a, Epsilon a) => SpMatrix a -> SpMatrix a -> SpVector a -> SpVector a 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 -- | used for Incomplete LU : remove entries in m corresponding to -- zero entries in m2 ilu0 :: (Epsilon a, Real a, Fractional a) => SpMatrix a -> (SpMatrix a, SpMatrix a) -- | `mSsor aa omega` : if `omega = 1` it returns the symmetric -- Gauss-Seidel preconditioner mSsor :: Fractional a => SpMatrix a -> a -> (SpMatrix a, SpMatrix a) -- | Partition a matrix into strictly subdiagonal, diagonal and strictly -- superdiagonal parts diagPartitions :: SpMatrix a -> (SpMatrix a, SpMatrix a, SpMatrix a) randArray :: PrimMonad m => Int -> Double -> Double -> m [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) -- | Sparsify an SpVector sparsifySV :: Epsilon a => SpVector a -> SpVector a -- | 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] -- | iterate until convergence is verified or we run out of a fixed -- iteration budget untilConverged :: MonadState a m => (a -> SpVector Double) -> (a -> a) -> m a diffSqL :: Floating a => [a] -> a 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