1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310
|
#!/usr/bin/env python
"""
PyCUDA-based FFT functions.
"""
import pycuda.driver as drv
import pycuda.gpuarray as gpuarray
import pycuda.elementwise as el
from pycuda.tools import context_dependent_memoize
import pycuda.tools as tools
import numpy as np
from . import cufft
from .cufft import CUFFT_COMPATIBILITY_NATIVE, \
CUFFT_COMPATIBILITY_FFTW_PADDING, \
CUFFT_COMPATIBILITY_FFTW_ASYMMETRIC, \
CUFFT_COMPATIBILITY_FFTW_ALL
from . import cudart
from . import misc
class Plan:
"""
CUFFT plan class.
This class represents an FFT plan for CUFFT.
Parameters
----------
shape : tuple of ints
Transform shape. May contain more than 3 elements.
in_dtype : { numpy.float32, numpy.float64, numpy.complex64, numpy.complex128 }
Type of input data.
out_dtype : { numpy.float32, numpy.float64, numpy.complex64, numpy.complex128 }
Type of output data.
batch : int
Number of FFTs to configure in parallel (default is 1).
stream : pycuda.driver.Stream
Stream with which to associate the plan. If no stream is specified,
the default stream is used.
mode : int
FFTW compatibility mode. Ignored in CUDA 9.2 and later.
inembed : numpy.array with dtype=numpy.int32
number of elements in each dimension of the input array
istride : int
distance between two successive input elements in the least significant
(innermost) dimension
idist : int
distance between the first element of two consective batches in the
input data
onembed : numpy.array with dtype=numpy.int32
number of elements in each dimension of the output array
ostride : int
distance between two successive output elements in the least significant
(innermost) dimension
odist : int
distance between the first element of two consective batches in the
output data
auto_allocate : bool
indicates whether the caller intends to allocate and manage the work area
"""
def __init__(self, shape, in_dtype, out_dtype, batch=1, stream=None,
mode=0x01, inembed=None, istride=1, idist=0, onembed=None,
ostride=1, odist=0, auto_allocate=True):
if np.isscalar(shape):
self.shape = (shape, )
else:
self.shape = shape
self.in_dtype = in_dtype
self.out_dtype = out_dtype
if batch <= 0:
raise ValueError('batch size must be greater than 0')
self.batch = batch
# Determine type of transformation:
if in_dtype == np.float32 and out_dtype == np.complex64:
self.fft_type = cufft.CUFFT_R2C
self.fft_func = cufft.cufftExecR2C
elif in_dtype == np.complex64 and out_dtype == np.float32:
self.fft_type = cufft.CUFFT_C2R
self.fft_func = cufft.cufftExecC2R
elif in_dtype == np.complex64 and out_dtype == np.complex64:
self.fft_type = cufft.CUFFT_C2C
self.fft_func = cufft.cufftExecC2C
elif in_dtype == np.float64 and out_dtype == np.complex128:
self.fft_type = cufft.CUFFT_D2Z
self.fft_func = cufft.cufftExecD2Z
elif in_dtype == np.complex128 and out_dtype == np.float64:
self.fft_type = cufft.CUFFT_Z2D
self.fft_func = cufft.cufftExecZ2D
elif in_dtype == np.complex128 and out_dtype == np.complex128:
self.fft_type = cufft.CUFFT_Z2Z
self.fft_func = cufft.cufftExecZ2Z
else:
raise ValueError('unsupported input/output type combination')
# Check for double precision support:
capability = misc.get_compute_capability(misc.get_current_device())
if capability < 1.3 and \
(misc.isdoubletype(in_dtype) or misc.isdoubletype(out_dtype)):
raise RuntimeError('double precision requires compute capability '
'>= 1.3 (you have %g)' % capability)
if inembed is not None:
inembed = inembed.ctypes.data
if onembed is not None:
onembed = onembed.ctypes.data
# Set up plan:
if len(self.shape) <= 0:
raise ValueError('invalid transform size')
n = np.asarray(self.shape, np.int32)
self.handle = cufft.cufftCreate()
# Set FFTW compatibility mode:
if cufft._cufft_version <= 9010:
cufft.cufftSetCompatibilityMode(self.handle, mode)
# Set auto-allocate mode
cufft.cufftSetAutoAllocation(self.handle, auto_allocate)
self.worksize = cufft.cufftMakePlanMany(
self.handle, len(self.shape), n.ctypes.data, inembed, istride, idist,
onembed, ostride, odist, self.fft_type, self.batch)
# Associate stream with plan:
if stream != None:
cufft.cufftSetStream(self.handle, stream.handle)
def set_work_area(self, work_area):
"""
Associate a caller-managed work area with the plan.
Parameters
----------
work_area : pycuda.gpuarray.GPUArray
"""
cufft.cufftSetWorkArea(self.handle, int(work_area.gpudata))
def __del__(self):
# Don't complain if handle destruction fails because the plan
# may have already been cleaned up:
try:
cufft.cufftDestroy(self.handle)
except:
pass
@context_dependent_memoize
def _get_scale_kernel(dtype):
ctype = tools.dtype_to_ctype(dtype)
return el.ElementwiseKernel(
"{ctype} scale, {ctype} *x".format(ctype=ctype),
"x[i] /= scale")
def _fft(x_gpu, y_gpu, plan, direction, scale=None):
"""
Fast Fourier Transform.
Parameters
----------
x_gpu : pycuda.gpuarray.GPUArray
Input array.
y_gpu : pycuda.gpuarray.GPUArray
Output array.
plan : Plan
FFT plan.
direction : { cufft.CUFFT_FORWARD, cufft.CUFFT_INVERSE }
Transform direction. Only affects in-place transforms.
Optional Parameters
-------------------
scale : int or float
Scale the values in the output array by dividing them by this value.
Notes
-----
This function should not be called directly.
"""
if (x_gpu.gpudata == y_gpu.gpudata) and \
plan.fft_type not in [cufft.CUFFT_C2C, cufft.CUFFT_Z2Z]:
raise ValueError('can only compute in-place transform of complex data')
if direction == cufft.CUFFT_FORWARD and \
plan.in_dtype in np.sctypes['complex'] and \
plan.out_dtype in np.sctypes['float']:
raise ValueError('cannot compute forward complex -> real transform')
if direction == cufft.CUFFT_INVERSE and \
plan.in_dtype in np.sctypes['float'] and \
plan.out_dtype in np.sctypes['complex']:
raise ValueError('cannot compute inverse real -> complex transform')
if plan.fft_type in [cufft.CUFFT_C2C, cufft.CUFFT_Z2Z]:
plan.fft_func(plan.handle, int(x_gpu.gpudata), int(y_gpu.gpudata),
direction)
else:
plan.fft_func(plan.handle, int(x_gpu.gpudata),
int(y_gpu.gpudata))
# Scale the result by dividing it by the number of elements:
if scale is not None:
func = _get_scale_kernel(y_gpu.dtype)
func(y_gpu.dtype.type(scale), y_gpu)
def fft(x_gpu, y_gpu, plan, scale=False):
"""
Fast Fourier Transform.
Compute the FFT of some data in device memory using the
specified plan.
Parameters
----------
x_gpu : pycuda.gpuarray.GPUArray
Input array.
y_gpu : pycuda.gpuarray.GPUArray
FFT of input array.
plan : Plan
FFT plan.
scale : bool, optional
If True, scale the computed FFT by the number of elements in
the input array.
Examples
--------
>>> import pycuda.autoinit
>>> import pycuda.gpuarray as gpuarray
>>> import numpy as np
>>> from skcuda.fft import fft, Plan
>>> N = 128
>>> x = np.asarray(np.random.rand(N), np.float32)
>>> xf = np.fft.fft(x)
>>> x_gpu = gpuarray.to_gpu(x)
>>> xf_gpu = gpuarray.empty(N/2+1, np.complex64)
>>> plan = Plan(x.shape, np.float32, np.complex64)
>>> fft(x_gpu, xf_gpu, plan)
>>> np.allclose(xf[0:N/2+1], xf_gpu.get(), atol=1e-6)
True
Returns
-------
y_gpu : pycuda.gpuarray.GPUArray
Computed FFT.
Notes
-----
For real to complex transformations, this function computes
N/2+1 non-redundant coefficients of a length-N input signal.
"""
if scale == True:
_fft(x_gpu, y_gpu, plan, cufft.CUFFT_FORWARD, x_gpu.size/plan.batch)
else:
_fft(x_gpu, y_gpu, plan, cufft.CUFFT_FORWARD)
def ifft(x_gpu, y_gpu, plan, scale=False):
"""
Inverse Fast Fourier Transform.
Compute the inverse FFT of some data in device memory using the
specified plan.
Parameters
----------
x_gpu : pycuda.gpuarray.GPUArray
Input array.
y_gpu : pycuda.gpuarray.GPUArray
Inverse FFT of input array.
plan : Plan
FFT plan.
scale : bool, optional
If True, scale the computed inverse FFT by the number of
elements in the output array.
Examples
--------
>>> import pycuda.autoinit
>>> import pycuda.gpuarray as gpuarray
>>> import numpy as np
>>> from skcuda.fft import fft, Plan
>>> N = 128
>>> x = np.asarray(np.random.rand(N), np.float32)
>>> xf = np.asarray(np.fft.fft(x), np.complex64)
>>> xf_gpu = gpuarray.to_gpu(xf[0:N/2+1])
>>> x_gpu = gpuarray.empty(N, np.float32)
>>> plan = Plan(N, np.complex64, np.float32)
>>> ifft(xf_gpu, x_gpu, plan, True)
>>> np.allclose(x, x_gpu.get(), atol=1e-6)
True
Notes
-----
For complex to real transformations, this function assumes the
input contains N/2+1 non-redundant FFT coefficents of a signal of
length N.
"""
if scale == True:
_fft(x_gpu, y_gpu, plan, cufft.CUFFT_INVERSE, y_gpu.size/plan.batch)
else:
_fft(x_gpu, y_gpu, plan, cufft.CUFFT_INVERSE)
if __name__ == "__main__":
import doctest
doctest.testmod()
|