A numpy.ndarray work-alike that stores its data and performs its computations on the compute device. shape and dtype work exactly as in numpy. Arithmetic methods in GPUArray support the broadcasting of scalars. (e.g. array+5) If the
allocator is a callable that, upon being called with an argument of the number of bytes to be allocated, returns an object that can be cast to an int representing the address of the newly allocated memory. Observe that both pycuda.driver.mem_alloc() and pycuda.tools.DeviceMemoryPool.alloc() are a model of this interface.
Returns the size of the leading dimension of self.
Transfer the contents the numpy.ndarray object ary onto the device.
ary must have the same dtype and size (not necessarily shape) as self.
Asynchronously transfer the contents the numpy.ndarray object ary onto the device, optionally sequenced on stream.
ary must have the same dtype and size (not necessarily shape) as self.
Return the real part of self, or self if it is real.
New in version 0.94.
Bind self to the pycuda.driver.TextureReference texref.
Due to alignment requirements, the effective texture bind address may be different from the requested one by an offset. This method returns this offset in units of self‘s data type. If allow_offset is False, a nonzero value of this offset will cause an exception to be raised.
Note
It is recommended to use bind_to_texref_ext() instead of this method.
Bind self to the pycuda.driver.TextureReference texref. In addition, set the texture reference’s format to match dtype and its channel count to channels. This routine also sets the texture reference’s pycuda.driver.TRSF_READ_AS_INTEGER flag, if necessary.
Due to alignment requirements, the effective texture bind address may be different from the requested one by an offset. This method returns this offset in units of self‘s data type. If allow_offset is False, a nonzero value of this offset will cause an exception to be raised.
New in version 0.93.
As of this writing, CUDA textures do not natively support double-precision floating point data. To remedy this deficiency, PyCUDA contains a workaround, which can be enabled by passing True for allow_double_hack. In this case, use the following code for texture access in your kernel code:
#include <pycuda-helpers.hpp>
texture<fp_tex_double, 1, cudaReadModeElementType> my_tex;
__global__ void f()
{
...
fp_tex1Dfetch(my_tex, threadIdx.x);
...
}
(This workaround was added in version 0.94.)
Return a GPUArray that is an exact copy of the numpy.ndarray instance ary.
See GPUArray for the meaning of allocator.
Return a GPUArray that is an exact copy of the numpy.ndarray instance ary. The copy is done asynchronously, optionally sequenced into stream.
See GPUArray for the meaning of allocator.
Create a GPUArray filled with numbers spaced step apart, starting from start and ending at stop.
For floating point arguments, the length of the result is ceil((stop - start)/step). This rule may result in the last element of the result being greater than stop.
dtype, if not specified, is taken as the largest common type of start, stop and step.
The pycuda.cumath module contains elementwise workalikes for the functions contained in math.
Warning
The following functionality is included in this documentation in the hope that it may be useful, but its interface may change in future revisions. Feedback is welcome.
Evaluating involved expressions on GPUArray instances can be somewhat inefficient, because a new temporary is created for each intermediate result. The functionality in the module pycuda.elementwise contains tools to help generate kernels that evaluate multi-stage expressions on one or several operands in a single pass.
Generate a kernel that takes a number of scalar or vector arguments and performs the scalar operation on each entry of its arguments, if that argument is a vector.
arguments is specified as a string formatted as a C argument list. operation is specified as a C assignment statement, without a semicolon. Vectors in operation should be indexed by the variable i.
name specifies the name as which the kernel is compiled, keep and options are passed unmodified to pycuda.compiler.SourceModule.
Here’s a usage example:
import pycuda.gpuarray as gpuarray
import pycuda.driver as cuda
import pycuda.autoinit
import numpy
from pycuda.curandom import rand as curand
a_gpu = curand((50,))
b_gpu = curand((50,))
from pycuda.elementwise import ElementwiseKernel
lin_comb = ElementwiseKernel(
"float a, float *x, float b, float *y, float *z",
"z[i] = a*x[i] + b*y[i]",
"linear_combination")
c_gpu = gpuarray.empty_like(a_gpu)
lin_comb(5, a_gpu, 6, b_gpu, c_gpu)
import numpy.linalg as la
assert la.norm((c_gpu - (5*a_gpu+6*b_gpu)).get()) < 1e-5
(You can find this example as examples/demo_elementwise.py in the PyCuda distribution.)
Generate a kernel that takes a number of scalar or vector arguments (at least one vector argument), performs the map_expr on each entry of the vector argument and then the reduce_expr on the outcome of that. neutral serves as an initial value. preamble offers the possibility to add preprocessor directives and other code (such as helper functions) to be added before the actual reduction kernel code.
Vectors in map_expr should be indexed by the variable i. reduce_expr uses the formal values “a” and “b” to indicate two operands of a binary reduction operation. If you do not specify a map_expr, “in[i]” – and therefore the presence of only one input argument – is automatically assumed.
dtype_out specifies the numpy.dtype in which the reduction is performed and in which the result is returned. neutral is specified as float or integer formatted as string. reduce_expr and map_expr are specified as string formatted operations and arguments is specified as a string formatted as a C argument list. name specifies the name as which the kernel is compiled, keep and options are passed unmodified to pycuda.compiler.SourceModule. preamble is specified as a string of code.
Here’s a usage example:
a = gpuarray.arange(400, dtype=numpy.float32)
b = gpuarray.arange(400, dtype=numpy.float32)
krnl = ReductionKernel(numpy.float32, neutral="0",
reduce_expr="a+b", map_expr="a[i]*b[i]",
arguments="float *a, float *b")
my_dot_prod = krnl(a, b).get()
Bogdan Opanchuk’s pycudafft offers a variety of GPU-based FFT implementations designed to work with pycuda.gpuarray.GPUArray objects.