Package fftw3f
[hide private]
[frames] | no frames]

Package fftw3f

source code


Python bindings to the FFTW library. 

Usage:

    In order to use PyFFTW you should have a basic understanding of the
    FFTW3 interface. For documentation about FFTW3 go to http://www.fftw.org
    In order to achieve maximum performance FFTW3 requires a planning stage
    where the actual FFT is created from an input and output array.
    To perform the FFT between the input and output array the plan is then
    executed. This interface is therefore significantly different from the
    traditional A = fft(B) interface.
    In contrast to the C-library PyFFTW utilizes the Plan class for planning.
    To create a fftw plan one creates a Plan object using an input and output
    array, and possible parameters. PyFFTW determines from the input and output
    arrays the correct plan to create. To perform the FFT one can either 
    call the Plan directly or call the method execute() or pass the plan
    to the execute function.
    

    Example:
    --------
    
    #create arrays
    >>>inputa = numpy.zeros((1024,3), dtype=complex)
    >>>outputa = numpy.zeros((1024,3), dtype=complex)
    
    # create a forward and backward fft plan
    >>>fft = fftw3.Plan(inputa,outputa, direction='forward', flags=['measure'])
    >>>ifft = fftw3.Plan(outputa, inputa, direction='backward', flags=['measure'])
    
    #initialize the input array
    >>>inputa[:] = 0
    >>>inputa += exp(-x**2/2)
    
    #perform a forward transformation
    >>>fft() # alternatively fft.execute() or fftw.execute(fft)
    
    # do some calculations with the output array
    >>>outputa *= D
    
    #perform a backward transformation
    >>>ifft() 
    
    The planning functions expect aligned, contiguous input arrays of 
    any shape. 
    Currently strides are not implemented. The dtype has to either be complex
    or double. If you want to perform ffts on single or longdouble precision
    arrays use the appropriate fftw3f or fftw3l module. FFTW overwrites the 
    arrays in the planning process, thus, if you use planning strategies 
    other than 'estimate' the arrays are going to be overwritten and have to 
    be reinitialized. 
    
    *IMPORTANT*
    -----------

    Because the plan uses pointers to the data of the arrays you cannot 
    perform operations on the arrays that change the data pointer. Therefore
    
    >>>a = zeros(1024, dtype=complex)
    >>>p = plan(a,b)
    >>>a = a+10
    >>>p()
    
    does not work, i.e. the a object references different memory, however 
    the Fourier transform will be performed on the original memory (the 
    plan actually contains a reference to the orgininal data (p.inarray), 
    otherwise this operation could even result in a python segfault).
    
    Aligned memory:
    ---------------
    
    On many platforms using the SIMD units for part of the floating point
    arithmetic significantly improves performance. FFTW can make use of the 
    SIMD operations, however the arrays have to be specifically aligned 
    in memory. PyFFTW provides a function which creates an numpy array 
    which is aligned to 
    a specified boundary. In most circumstances the default alignment to 16 
    byte boundary is what you want. Note that the same precautions as above     
    apply, i.e. creating an aligned array and then doing something like 
    a=a+1 will result in new memory allocated by python which might not 
    be aligned.
    
    PyFFTW interface naming conventions:
    ------------------------------------

    All exposed fftw-functions do have the same names as the C-functions with the
    leading fftw_ striped from the name.
    Direct access to the C-functions is available by importing lib.lib, the usual
    precautions for using C-functions from Python apply. 
    
    Advanced and Guru interface:
    ----------------------------

    Currently only the execute_dft function from the fftw guru and advanced 
    interface is exposed.
    It is explicitly name guru_execute_dft. You should only use
    these if you know what you're doing, as no checking is done on these functions.
    


Constants:

    fftw_flags      -- dictionary of possible flags for creating plans
    fft_direction   -- the direction of the fft (see the fftw documentation
                       for the mathematical meaning).
    realfft_type    -- a dictionary of possible types for real-to-real
                       transforms (see the fftw documentation for a 
                       more detailed description).

Submodules [hide private]

Classes [hide private]
  Plan
Object representing a fftw plan used to execute Fourier transforms in fftw
Functions [hide private]
 
export_wisdom_to_file(filename)
Export accumulated wisdom to file given by the filename
source code
 
export_wisdom_to_string()
Returns a string with the accumulated wisdom
source code
 
import_wisdom_from_string(wisdom)
Import wisdom from the given string
source code
 
import_wisdom_from_file(filename)
Imports wisdom from the file given by the filename
source code
 
import_system_wisdom()
Import the system wisdom, this lives under /etc/fftw/wisdom on Unix/Linux systems
source code
 
forget_wisdom()
Clear all wisdom
source code
 
create_aligned_array(shape, dtype=<type 'numpy.complex64'>, boundary=16)
Create an array which is aligned in memory
source code
 
execute(plan)
Execute fftw-plan, i.e.
source code
 
guru_execute_dft(plan, inarray, outarray)
Guru interface: perform Fourier transform on two arrays, outarray=fft(inarray) using the given plan.
source code
 
destroy_plan(plan)
Delete the given plan
source code
Variables [hide private]
  fftw_flags = {'conserve memory': 4, 'destroy input': 1, 'estim...
  fft_direction = {'backward': 1, 'forward': -1}
  realfft_type = {'discrete hartley': 2, 'halfcomplex c2r': 1, '...
  __package__ = 'fftw3f'
Function Details [hide private]

create_aligned_array(shape, dtype=<type 'numpy.complex64'>, boundary=16)

source code 
Create an array which is aligned in memory

Parameters:
    shape       --  the shape of the array 
    dtype       --  the dtype of the array (default=typeDict['singlecomplex'])
    boundary    --  the byte boundary to align to (default=16)

execute(plan)

source code 

Execute fftw-plan, i.e. perform Fourier transform on the arrays given when the plan was created

guru_execute_dft(plan, inarray, outarray)

source code 

Guru interface: perform Fourier transform on two arrays, outarray=fft(inarray) using the given plan. Important: This function does not perform any checks on the array shape and alignment for performance reasons. It is therefore crucial to only provide arrays with the same shape, dtype and alignment as the arrays used for planning, failure to do so can lead to unexpected behaviour and even python segfaulting.


Variables Details [hide private]

fftw_flags

Value:
{'conserve memory': 4,
 'destroy input': 1,
 'estimate': 64,
 'exhaustive': 8,
 'measure': 0,
 'patient': 32,
 'preserve input': 16,
 'unaligned': 2}

realfft_type

Value:
{'discrete hartley': 2,
 'halfcomplex c2r': 1,
 'halfcomplex r2c': 0,
 'realeven 00': 3,
 'realeven 01': 4,
 'realeven 10': 5,
 'realeven 11': 6,
 'realodd 00': 7,
...