fft-0.1.8.5: Bindings to the FFTW library.

Safe HaskellNone
LanguageHaskell98

Math.FFT.Base

Synopsis

Documentation

class (Storable a, RealFloat a) => FFTWReal a where Source #

Our API is polymorphic over the real data type. FFTW, at least in principle, supports single precision Float, double precision Double and long double CLDouble (presumable?).

lock :: MVar () Source #

This lock must be taken during planning of any transform. The FFTW library is not thread-safe in the planning phase. Thankfully, the lock is not needed during the execute phase.

withLock :: IO a -> IO a Source #

newtype Flag Source #

The Flag type is used to influence the kind of plans which are created. To specify multiple flags, use a bitwise .|..

Constructors

Flag 

Fields

Instances

nullFlag :: Flag Source #

Default flag. For most transforms, this is equivalent to setting measure and preserveInput. The exceptions are complex to real and half-complex to real transforms.

destroyInput :: Flag Source #

Allows FFTW to overwrite the input array with arbitrary data; this can sometimes allow more efficient algorithms to be employed.

Setting this flag implies that two memory allocations will be done, one for work space, and one for the result. When estimate is not set, we will be doing two memory allocations anyway, so we set this flag as well (since we don't retain the work array anyway).

preserveInput :: Flag Source #

preserveInput specifies that an out-of-place transform must not change its input array. This is ordinarily the default, except for complex to real transforms for which destroyInput is the default. In the latter cases, passing preserveInput will attempt to use algorithms that do not destroy the input, at the expense of worse performance; for multi-dimensional complex to real transforms, however, no input-preserving algorithms are implemented so the Haskell bindings will set destroyInput and do a transform with two memory allocations.

unaligned :: Flag Source #

Instruct FFTW not to generate a plan which uses SIMD instructions, even if the memory you are planning with is aligned. This should only be needed if you are using the guru interface and want to reuse a plan with memory that may be unaligned (i.e. you constructed the CArray with unsafeForeignPtrToCArray).

conserveMemory :: Flag Source #

The header claims that this flag is documented, but in reality, it is not. I don't know what it does and it is here only for completeness.

estimate :: Flag Source #

estimate specifies that, instead of actual measurements of different algorithms, a simple heuristic is used to pick a (probably sub-optimal) plan quickly. With this flag, the input/output arrays are not overwritten during planning.

This is the only planner flag for which a single memory allocation is possible.

measure :: Flag Source #

measure tells FFTW to find an optimized plan by actually computing several FFTs and measuring their execution time. Depending on your machine, this can take some time (often a few seconds). measure is the default planning option.

patient :: Flag Source #

patient is like measure, but considers a wider range of algorithms and often produces a "more optimal" plan (especially for large transforms), but at the expense of several times longer planning time (especially for large transforms).

exhaustive :: Flag Source #

exhaustive is like patient but considers an even wider range of algorithms, including many that we think are unlikely to be fast, to produce the most optimal plan but with a substantially increased planning time.

data Sign Source #

Determine which direction of DFT to execute.

Constructors

DFTForward 
DFTBackward 

Instances

Eq Sign Source # 

Methods

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

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

Show Sign Source # 

Methods

showsPrec :: Int -> Sign -> ShowS #

show :: Sign -> String #

showList :: [Sign] -> ShowS #

data Kind Source #

Real to Real transform kinds.

Instances

Eq Kind Source # 

Methods

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

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

Show Kind Source # 

Methods

showsPrec :: Int -> Kind -> ShowS #

show :: Kind -> String #

showList :: [Kind] -> ShowS #

type TSpec = ([IODim], [IODim]) Source #

Tuple of transform dimensions and non-transform dimensions of the array.

data DFT Source #

Types of transforms. Used to control dftShape.

Constructors

CC 
RC 
CR 
CRO 
RR 

Instances

Eq DFT Source # 

Methods

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

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

Show DFT Source # 

Methods

showsPrec :: Int -> DFT -> ShowS #

show :: DFT -> String #

showList :: [DFT] -> ShowS #

check :: Plan -> IO () Source #

Verify that a plan is valid. Thows an exception if not.

execute :: Plan -> IO () Source #

Confirm that the plan is valid, then execute the transform.

unsafeNormalize :: (Ix i, Shapable i, Fractional e, Storable e) => [Int] -> CArray i e -> CArray i e Source #

In-place normalization outside of IO. You must be able to prove that no reference to the original can be retained.

dftG :: (FFTWReal r, Ix i, Shapable i) => Sign -> Flag -> [Int] -> CArray i (Complex r) -> CArray i (Complex r) Source #

Normalized general complex DFT

dftCRG :: (FFTWReal r, Ix i, Shapable i) => Flag -> [Int] -> CArray i (Complex r) -> CArray i r Source #

Normalized general complex to real DFT where the last transformed dimension is logically even.

dftCROG :: (FFTWReal r, Ix i, Shapable i) => Flag -> [Int] -> CArray i (Complex r) -> CArray i r Source #

Normalized general complex to real DFT where the last transformed dimension is logicall odd.

dftN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i (Complex r) -> CArray i (Complex r) Source #

Multi-dimensional forward DFT.

idftN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i (Complex r) -> CArray i (Complex r) Source #

Multi-dimensional inverse DFT.

dftRCN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i (Complex r) Source #

Multi-dimensional forward DFT of real data.

dftCRN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i (Complex r) -> CArray i r Source #

Multi-dimensional inverse DFT of Hermitian-symmetric data (where only the non-negative frequencies are given).

dftCRON :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i (Complex r) -> CArray i r Source #

Multi-dimensional inverse DFT of Hermitian-symmetric data (where only the non-negative frequencies are given) and the last transformed dimension is logically odd.

fzr :: b -> [a] -> [(a, b)] Source #

drr :: (FFTWReal r, Ix i, Shapable i) => Kind -> [Int] -> CArray i r -> CArray i r Source #

dftRRN :: (FFTWReal r, Ix i, Shapable i) => [(Int, Kind)] -> CArray i r -> CArray i r Source #

Multi-dimensional real to real transform. The result is not normalized.

dftRHN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional real to half-complex transform. The result is not normalized.

dftHRN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional half-complex to real transform. The result is not normalized.

dhtN :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Discrete Hartley Transform. The result is not normalized.

dct1N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 1 discrete cosine transform.

dct2N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 2 discrete cosine transform. This is commonly known as the DCT.

dct3N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 3 discrete cosine transform. This is commonly known as the inverse DCT. The result is not normalized.

dct4N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 4 discrete cosine transform.

dst1N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 1 discrete sine transform.

dst2N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 2 discrete sine transform.

dst3N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 3 discrete sine transform.

dst4N :: (FFTWReal r, Ix i, Shapable i) => [Int] -> CArray i r -> CArray i r Source #

Multi-dimensional Type 4 discrete sine transform.

dft :: (FFTWReal r, Ix i, Shapable i) => CArray i (Complex r) -> CArray i (Complex r) Source #

1-dimensional complex DFT.

idft :: (FFTWReal r, Ix i, Shapable i) => CArray i (Complex r) -> CArray i (Complex r) Source #

1-dimensional complex inverse DFT. Inverse of dft.

dftRC :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i (Complex r) Source #

1-dimensional real to complex DFT.

dftCR :: (FFTWReal r, Ix i, Shapable i) => CArray i (Complex r) -> CArray i r Source #

1-dimensional complex to real DFT with logically even dimension. Inverse of dftRC.

dftCRO :: (FFTWReal r, Ix i, Shapable i) => CArray i (Complex r) -> CArray i r Source #

1-dimensional complex to real DFT with logically odd dimension. Inverse of dftRC.

dftRH :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional real to half-complex DFT.

dftHR :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional half-complex to real DFT. Inverse of dftRH after normalization.

dht :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Discrete Hartley Transform. Self-inverse after normalization.

dct1 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 1 discrete cosine transform.

dct2 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 2 discrete cosine transform. This is commonly known as the DCT.

dct3 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 3 discrete cosine transform. This is commonly known as the inverse DCT.

dct4 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 4 discrete cosine transform.

dst1 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 1 discrete sine transform.

dst2 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 2 discrete sine transform.

dst3 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 3 discrete sine transform.

dst4 :: (FFTWReal r, Ix i, Shapable i) => CArray i r -> CArray i r Source #

1-dimensional Type 4 discrete sine transform.

has :: Flag -> Flag -> Bool infix 7 Source #

transformCArray :: (Ix i, Storable a, Storable b) => Flag -> CArray i a -> (i, i) -> (FFTWFlag -> Ptr a -> Ptr b -> IO Plan) -> CArray i b Source #

Try to transform a CArray with only one memory allocation (for the result). If we can find a way to prove that FFTW already has a sufficiently good plan for this transform size and the input will not be overwritten, then we could call have a version of this that does not require estimate. Since this is not currently the case, we require estimate to be set. Note that we do not check for the preserveInput flag here. This is because the default is to preserve input for all but the C->R and HC->R transforms. Therefore, this function must not be called for those transforms, unless preserveInput is set.

transformCArray' :: (Ix i, Storable a, Storable b) => Flag -> CArray i a -> (i, i) -> (FFTWFlag -> Ptr a -> Ptr b -> IO Plan) -> CArray i b Source #

Transform a CArray with two memory allocations. This is entirely safe with all transforms, but it must allocate a temporary array to do the planning in.

dftShape :: (Ix i, Shapable i, Storable e) => DFT -> [Int] -> CArray i e -> ((i, i), TSpec) Source #

All the logic for determining shape of resulting array, and how to do the transform.

withTSpec :: TSpec -> (CInt -> Ptr IODim -> CInt -> Ptr IODim -> IO a) -> IO a Source #

A simple helper.

adjust :: (a -> a) -> Int -> [a] -> [a] Source #

A generally useful list utility

dftGU :: (FFTWReal r, Ix i, Shapable i) => Sign -> Flag -> [Int] -> CArray i (Complex r) -> CArray i (Complex r) Source #

Complex to Complex DFT, un-normalized.

dftRCG :: (FFTWReal r, Ix i, Shapable i) => Flag -> [Int] -> CArray i r -> CArray i (Complex r) Source #

Real to Complex DFT.

dftCRG_ :: (FFTWReal r, Ix i, Shapable i) => Bool -> Flag -> [Int] -> CArray i (Complex r) -> CArray i r Source #

Complex to Real DFT. The first argument determines whether the last transformed dimension is logically odd or even. True implies the dimension is odd.

dftCRGU :: (FFTWReal r, Ix i, Shapable i) => Flag -> [Int] -> CArray i (Complex r) -> CArray i r Source #

Complex to Real DFT where last transformed dimension is logically even.

dftCROGU :: (FFTWReal r, Ix i, Shapable i) => Flag -> [Int] -> CArray i (Complex r) -> CArray i r Source #

Complex to Real DFT where last transformed dimension is logically odd.

dftRRG :: (FFTWReal r, Ix i, Shapable i) => Flag -> [(Int, Kind)] -> CArray i r -> CArray i r Source #

Real to Real transforms.

exportWisdomString :: IO String Source #

Queries the FFTW cache. The String can be written to a file so the wisdom can be reused on a subsequent run.

importWisdomString :: String -> IO Bool Source #

Add wisdom to the FFTW cache. Returns True if it is successful.

importWisdomSystem :: IO Bool Source #

Tries to import wisdom from a global source, typically etcfftw/wisdom. Returns True if it was successful.

c_rodft10 :: FFTWKind Source #

Corresponds to the fftw_iodim structure. It completely describes the layout of each dimension, before and after the transform.

data IODim Source #

Constructors

IODim 

Fields

Instances

Eq IODim Source # 

Methods

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

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

Data IODim Source # 

Methods

gfoldl :: (forall d b. Data d => c (d -> b) -> d -> c b) -> (forall g. g -> c g) -> IODim -> c IODim #

gunfold :: (forall b r. Data b => c (b -> r) -> c r) -> (forall r. r -> c r) -> Constr -> c IODim #

toConstr :: IODim -> Constr #

dataTypeOf :: IODim -> DataType #

dataCast1 :: Typeable (* -> *) t => (forall d. Data d => c (t d)) -> Maybe (c IODim) #

dataCast2 :: Typeable (* -> * -> *) t => (forall d e. (Data d, Data e) => c (t d e)) -> Maybe (c IODim) #

gmapT :: (forall b. Data b => b -> b) -> IODim -> IODim #

gmapQl :: (r -> r' -> r) -> r -> (forall d. Data d => d -> r') -> IODim -> r #

gmapQr :: (r' -> r -> r) -> r -> (forall d. Data d => d -> r') -> IODim -> r #

gmapQ :: (forall d. Data d => d -> u) -> IODim -> [u] #

gmapQi :: Int -> (forall d. Data d => d -> u) -> IODim -> u #

gmapM :: Monad m => (forall d. Data d => d -> m d) -> IODim -> m IODim #

gmapMp :: MonadPlus m => (forall d. Data d => d -> m d) -> IODim -> m IODim #

gmapMo :: MonadPlus m => (forall d. Data d => d -> m d) -> IODim -> m IODim #

Show IODim Source # 

Methods

showsPrec :: Int -> IODim -> ShowS #

show :: IODim -> String #

showList :: [IODim] -> ShowS #

Storable IODim Source # 

Methods

sizeOf :: IODim -> Int #

alignment :: IODim -> Int #

peekElemOff :: Ptr IODim -> Int -> IO IODim #

pokeElemOff :: Ptr IODim -> Int -> IODim -> IO () #

peekByteOff :: Ptr b -> Int -> IO IODim #

pokeByteOff :: Ptr b -> Int -> IODim -> IO () #

peek :: Ptr IODim -> IO IODim #

poke :: Ptr IODim -> IODim -> IO () #

type Plan = Ptr FFTWPlan Source #

A plan is an opaque foreign object.

type FFTWPlan = () Source #

cf_plan_guru_dft :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr (Complex Float) -> Ptr (Complex Float) -> FFTWSign -> FFTWFlag -> IO Plan Source #

Plan a complex to complex transform using the guru interface.

cf_plan_guru_dft_r2c :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr Float -> Ptr (Complex Float) -> FFTWFlag -> IO Plan Source #

Plan a real to complex transform using the guru interface.

cf_plan_guru_dft_c2r :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr (Complex Float) -> Ptr Float -> FFTWFlag -> IO Plan Source #

Plan a complex to real transform using the guru interface.

cf_plan_guru_r2r :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr Float -> Ptr Float -> Ptr FFTWKind -> FFTWFlag -> IO Plan Source #

Plan a real to real transform using the guru interface.

c_plan_guru_dft :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr (Complex Double) -> Ptr (Complex Double) -> FFTWSign -> FFTWFlag -> IO Plan Source #

Plan a complex to complex transform using the guru interface.

c_plan_guru_dft_r2c :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr Double -> Ptr (Complex Double) -> FFTWFlag -> IO Plan Source #

Plan a real to complex transform using the guru interface.

c_plan_guru_dft_c2r :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr (Complex Double) -> Ptr Double -> FFTWFlag -> IO Plan Source #

Plan a complex to real transform using the guru interface.

c_plan_guru_r2r :: CInt -> Ptr IODim -> CInt -> Ptr IODim -> Ptr Double -> Ptr Double -> Ptr FFTWKind -> FFTWFlag -> IO Plan Source #

Plan a real to real transform using the guru interface.

c_execute :: Plan -> IO () Source #

Simple plan execution

c_free :: Ptr a -> IO () Source #

Frees memory allocated by fftw_malloc. Currently, we only need this to free the wisdom string.