module Numerical.HBLAS.BLAS(
GemvFun
,GemmFun
,TrsvFun
,dgemm
,sgemm
,cgemm
,zgemm
,sgemv
,dgemv
,cgemv
,zgemv
,strsv
,dtrsv
,ctrsv
,ztrsv
) where
import Numerical.HBLAS.UtilsFFI
import Numerical.HBLAS.BLAS.FFI
import Numerical.HBLAS.BLAS.Internal
import Numerical.HBLAS.MatrixTypes
import Control.Monad.Primitive
import Data.Complex
type GemmFun el orient s m = Transpose ->Transpose -> el -> el -> MDenseMatrix s orient el
-> MDenseMatrix s orient el -> MDenseMatrix s orient el -> m ()
type GemvFun el orient s m = Transpose -> el -> el
-> MDenseMatrix s orient el -> MDenseVector s Direct el -> MDenseVector s Direct el -> m ()
type TrsvFun el orient s m =
MatUpLo -> Transpose -> MatDiag
-> MDenseMatrix s orient el -> MDenseVector s Direct el -> m ()
sgemm :: PrimMonad m=> GemmFun Float orient (PrimState m) m
sgemm = gemmAbstraction "sgemm" cblas_sgemm_safe cblas_sgemm_unsafe (\x f -> f x )
dgemm :: PrimMonad m=> GemmFun Double orient (PrimState m) m
dgemm = gemmAbstraction "dgemm" cblas_dgemm_safe cblas_dgemm_unsafe (\x f -> f x )
cgemm :: PrimMonad m=> GemmFun (Complex Float) orient (PrimState m) m
cgemm = gemmAbstraction "cgemm" cblas_cgemm_safe cblas_cgemm_unsafe withRStorable_
zgemm :: PrimMonad m=> GemmFun (Complex Double) orient (PrimState m) m
zgemm = gemmAbstraction "zgemm" cblas_zgemm_safe cblas_zgemm_unsafe withRStorable_
sgemv :: PrimMonad m => GemvFun Float orient (PrimState m) m
sgemv = gemvAbstraction "sgemv" cblas_sgemv_safe cblas_sgemv_unsafe (flip ($))
dgemv :: PrimMonad m => GemvFun Double orient (PrimState m) m
dgemv = gemvAbstraction "dgemv" cblas_dgemv_safe cblas_dgemv_unsafe (flip ($))
cgemv :: PrimMonad m => GemvFun (Complex Float) orient (PrimState m) m
cgemv = gemvAbstraction "cgemv" cblas_cgemv_safe cblas_cgemv_unsafe withRStorable_
zgemv :: PrimMonad m => GemvFun (Complex Double) orient (PrimState m) m
zgemv = gemvAbstraction "zgemv" cblas_zgemv_safe cblas_zgemv_unsafe withRStorable_
strsv :: PrimMonad m => TrsvFun Float orient (PrimState m) m
strsv = trsvAbstraction "strsv" cblas_strsv_safe cblas_strsv_unsafe
dtrsv :: PrimMonad m => TrsvFun Double orient (PrimState m) m
dtrsv = trsvAbstraction "dtrsv" cblas_dtrsv_safe cblas_dtrsv_unsafe
ctrsv :: PrimMonad m => TrsvFun (Complex Float) orient (PrimState m) m
ctrsv = trsvAbstraction "ctrsv" cblas_ctrsv_safe cblas_ctrsv_unsafe
ztrsv :: PrimMonad m => TrsvFun (Complex Double) orient (PrimState m) m
ztrsv = trsvAbstraction "ztrsv" cblas_ztrsv_safe cblas_ztrsv_unsafe