GPU Arrays#

Vector Types#

class pycuda.gpuarray.vec#

All of CUDA’s supported vector types, such as float3 and long4 are available as numpy data types within this class. These numpy.dtype instances have field names of x, y, z, and w just like their CUDA counterparts. They will work both for parameter passing to kernels as well as for passing data back and forth between kernels and Python code. For each type, a make_type function is also provided (e.g. make_float3(x,y,z)).

The GPUArray Array Class#

class pycuda.gpuarray.GPUArray(shape, dtype, *, allocator=None, order='C')#

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.

All arguments beyond allocator should be considered keyword-only.

gpudata#

The pycuda.driver.DeviceAllocation instance created for the memory that backs this GPUArray.

shape#

The tuple of lengths of each dimension in the array.

dtype#

The numpy.dtype of the items in the GPU array.

size#

The number of meaningful entries in the array. Can also be computed by multiplying up the numbers in shape.

mem_size#

The total number of entries, including padding, that are present in the array. Padding may arise for example because of pitch adjustment by pycuda.driver.mem_alloc_pitch().

nbytes#

The size of the entire array in bytes. Computed as size times dtype.itemsize.

strides#

Tuple of bytes to step in each dimension when traversing an array.

flags#

Return an object with attributes c_contiguous, f_contiguous and forc, which may be used to query contiguity properties in analogy to numpy.ndarray.flags.

ptr#

Return an int reflecting the address in device memory where this array resides.

__cuda_array_interface__#

Return a CUDA Array Interface dict describing this array’s data.

__len__()#

Returns the size of the leading dimension of self.

Warning

This method existed in version 0.93 and below, but it returned the value of size instead of its current value. The change was made in order to match numpy.

reshape(shape, order='C')#

Returns an array containing the same data with a new shape.

ravel()#

Returns flattened array containing the same data.

view(dtype=None)#

Returns view of array with the same data. If dtype is different from current dtype, the actual bytes of memory will be reinterpreted.

squeeze(dtype=None)#

Returns a view of the array with dimensions of length 1 removed.

set(ary)#

Transfer the contents the numpy.ndarray object ary onto the device.

ary must have the same dtype and size (not necessarily shape) as self.

set_async(ary, stream=None)#

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.

get(ary=None, pagelocked=False)#

Transfer the contents of self into ary or a newly allocated numpy.ndarray. If ary is given, it must have the same shape and dtype. If it is not given, a pagelocked specifies whether the new array is allocated page-locked.

Changed in version 2015.2: ary with different shape was deprecated.

get_async(stream=None, ary=None)#

Transfer the contents of self into ary or a newly allocated numpy.ndarray. If ary is given, it must have the right size (not necessarily shape) and dtype. If it is not given, a page-locked array is newly allocated.

copy()#

New in version 2013.1.

mul_add(self, selffac, other, otherfac, add_timer=None, stream=None):

Return selffac*self + otherfac*other. add_timer, if given, is invoked with the result from pycuda.driver.Function.prepared_timed_call().

__add__(other)#
__sub__(other)#
__iadd__(other)#
__isub__(other)#
__neg__(other)#
__mul__(other)#
__div__(other)#
__rdiv__(other)#
__pow__(other)#
__abs__()#

Return a GPUArray containing the absolute value of each element of self.

fill(scalar, stream=None)#

Fill the array with scalar.

astype(dtype, stream=None)#

Return self, cast to dtype.

any(stream=None, allocator=None)#
all(stream=None, allocator=None)#
real#

Return the real part of self, or self if it is real.

New in version 0.94.

imag#

Return the imaginary part of self, or zeros_like(self) if it is real.

conj(out=None)#

Return the complex conjugate of self, or self if it is real. If out is not given, a newly allocated GPUArray will returned. Use out=self to get conjugate in-place.

Changed in version 2020.1.1: add out parameter

conjugate(out=None)#

alias of conj()

New in version 2020.1.1.

bind_to_texref(texref, allow_offset=False)#

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_to_texref_ext(texref, channels=1, allow_double_hack=False, allow_offset=False)#

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.)

Constructing GPUArray Instances#

pycuda.gpuarray.to_gpu(ary, allocator=None)#

Return a GPUArray that is an exact copy of the numpy.ndarray instance ary.

See GPUArray for the meaning of allocator.

pycuda.gpuarray.to_gpu_async(ary, allocator=None, stream=None)#

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.

pycuda.gpuarray.empty(shape, dtype, *, allocator=None, order='C')#

A synonym for the GPUArray constructor.

pycuda.gpuarray.zeros(shape, dtype=np.float64, *, allocator=None, order='C')#

Same as empty(), but the GPUArray is zero-initialized before being returned.

pycuda.gpuarray.ones(shape, dtype=np.float64, *, allocator=None, order='C')#

Same as empty(), but the GPUArray is one-initialized before being returned.

pycuda.gpuarray.empty_like(other_ary, dtype=None, order='K')#

Make a new, uninitialized GPUArray having the same properties as other_ary. The dtype and order attributes allow these aspects to be set independently of their values in other_ary. For order, “A” means retain Fortran-ordering if the input is Fortran-contiguous, otherwise use “C” ordering. The default, order or “K” tries to match the strides of other_ary as closely as possible.

pycuda.gpuarray.zeros_like(other_ary, dtype=None, order='K')#

Make a new, zero-initialized GPUArray having the same properties as other_ary. The dtype and order attributes allow these aspects to be set independently of their values in other_ary. For order, “A” means retain Fortran-ordering if the input is Fortran-contiguous, otherwise use “C” ordering. The default, order or “K” tries to match the strides of other_ary as closely as possible.

pycuda.gpuarray.ones_like(other_ary, dtype=None, order='K')#

Make a new, ones-initialized GPUArray having the same properties as other_ary. The dtype and order attributes allow these aspects to be set independently of their values in other_ary. For order, “A” means retain Fortran-ordering if the input is Fortran-contiguous, otherwise use “C” ordering. The default, order or “K” tries to match the strides of other_ary as closely as possible.

pycuda.gpuarray.arange(start, stop, step, dtype=None, stream=None)#

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.

pycuda.gpuarray.take(a, indices, stream=None)#

Return the GPUArray [a[indices[0]], ..., a[indices[n]]]. For the moment, a must be a type that can be bound to a texture.

pycuda.gpuarray.concatenate(arrays, axis=0, allocator=None)#

Join a sequence of arrays along an existing axis.

pycuda.gpuarray.stack(arrays, axis=0, allocator=None)#

Join a sequence of arrays along a new axis.

pycuda.gpuarray.logical_and(x1, x2, /, out=None, * allocator=None)#

Returns the elementwise logical AND values of x1 and x2.

pycuda.gpuarray.logical_or(x1, x2, /, out=None, * allocator=None)#

Returns the elementwise logical OR values of x1 and x2.

pycuda.gpuarray.logical_not(x, /, out=None, * allocator=None)#

Returns the elementwise logical NOT of x.

Conditionals#

pycuda.gpuarray.if_positive(criterion, then_, else_, out=None, stream=None)#

Return an array like then_, which, for the element at index i, contains then_[i] if criterion[i]>0, else else_[i]. (added in 0.94)

pycuda.gpuarray.maximum(a, b, out=None, stream=None)#

Return the elementwise maximum of a and b. (added in 0.94)

pycuda.gpuarray.minimum(a, b, out=None, stream=None)#

Return the elementwise minimum of a and b. (added in 0.94)

Reductions#

pycuda.gpuarray.sum(a, dtype=None, stream=None)#
pycuda.gpuarray.any(a, stream=None, allocator=None)#
pycuda.gpuarray.all(a, stream=None, allocator=None)#
pycuda.gpuarray.subset_sum(subset, a, dtype=None, stream=None)#

New in version 2013.1.

pycuda.gpuarray.dot(a, b, dtype=None, stream=None)#
pycuda.gpuarray.subset_dot(subset, a, b, dtype=None, stream=None)#
pycuda.gpuarray.max(a, stream=None)#
pycuda.gpuarray.min(a, stream=None)#
pycuda.gpuarray.subset_max(subset, a, stream=None)#
pycuda.gpuarray.subset_min(subset, a, stream=None)#

Elementwise Functions on GPUArray Instances#

The pycuda.cumath module contains elementwise workalikes for the functions contained in math.

Rounding and Absolute Value#

pycuda.cumath.fabs(array, *, out=None, stream=None)#
pycuda.cumath.ceil(array, *, out=None, stream=None)#
pycuda.cumath.floor(array, *, out=None, stream=None)#

Exponentials, Logarithms and Roots#

pycuda.cumath.exp(array, *, out=None, stream=None)#
pycuda.cumath.log(array, *, out=None, stream=None)#
pycuda.cumath.log10(array, *, out=None, stream=None)#
pycuda.cumath.sqrt(array, *, out=None, stream=None)#

Trigonometric Functions#

pycuda.cumath.sin(array, *, out=None, stream=None)#
pycuda.cumath.cos(array, *, out=None, stream=None)#
pycuda.cumath.tan(array, *, out=None, stream=None)#
pycuda.cumath.asin(array, *, out=None, stream=None)#
pycuda.cumath.acos(array, *, out=None, stream=None)#
pycuda.cumath.atan(array, *, out=None, stream=None)#

Hyperbolic Functions#

pycuda.cumath.sinh(array, *, out=None, stream=None)#
pycuda.cumath.cosh(array, *, out=None, stream=None)#
pycuda.cumath.tanh(array, *, out=None, stream=None)#

Floating Point Decomposition and Assembly#

pycuda.cumath.fmod(arg, mod, stream=None)#

Return the floating point remainder of the division arg/mod, for each element in arg and mod.

pycuda.cumath.frexp(arg, stream=None)#

Return a tuple (significands, exponents) such that arg == significand * 2**exponent.

pycuda.cumath.ldexp(significand, exponent, stream=None)#

Return a new array of floating point values composed from the entries of significand and exponent, paired together as result = significand * 2**exponent.

pycuda.cumath.modf(arg, stream=None)#

Return a tuple (fracpart, intpart) of arrays containing the integer and fractional parts of arg.

Generating Arrays of Random Numbers#

pycuda.curandom.rand(shape, dtype=numpy.float32, stream=None)#

Return an array of shape filled with random values of dtype in the range [0,1).

Note

The use case for this function is “I need some random numbers. I don’t care how good they are or how fast I get them.” It uses a pretty terrible MD5-based generator and doesn’t even attempt to cache generated code.

If you’re interested in a non-toy random number generator, use the CURAND-based functionality below.

Warning

The following classes are using random number generators that run on the GPU. Each thread uses its own generator. Creation of those generators requires more resources than subsequent generation of random numbers. After experiments it looks like maximum number of active generators on Tesla devices (with compute capabilities 1.x) is 256. Fermi devices allow for creating 1024 generators without any problems. If there are troubles with creating objects of class PseudoRandomNumberGenerator or QuasiRandomNumberGenerator decrease number of created generators (and therefore number of active threads).

A pseudorandom sequence of numbers satisfies most of the statistical properties of a truly random sequence but is generated by a deterministic algorithm. A quasirandom sequence of n-dimensional points is generated by a deterministic algorithm designed to fill an n-dimensional space evenly.

Quasirandom numbers are more expensive to generate.

pycuda.curandom.get_curand_version()#

Obtain the version of CURAND against which PyCUDA was compiled. Returns a 3-tuple of integers as (major, minor, revision).

pycuda.curandom.seed_getter_uniform(N)#

Return an GPUArray filled with one random int32 repeated N times which can be used as a seed for XORWOW generator.

pycuda.curandom.seed_getter_unique(N)#

Return an GPUArray filled with N random int32 which can be used as a seed for XORWOW generator.

class pycuda.curandom.XORWOWRandomNumberGenerator(seed_getter=None, offset=0)#
Parameters:
  • seed_getter – a function that, given an integer count, will yield an int32 GPUArray of seeds.

  • offset – Starting index into the XORWOW sequence, given seed.

Provides pseudorandom numbers. Generates sequences with period at least \(2^190\).

CUDA 3.2 and above.

New in version 2011.1.

fill_uniform(data, stream=None)#

Fills in GPUArray data with uniformly distributed pseudorandom values.

gen_uniform(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with uniformly distributed pseudorandom values, and returns newly created object.

fill_normal(data, stream=None)#

Fills in GPUArray data with normally distributed pseudorandom values.

gen_normal(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with normally distributed pseudorandom values, and returns newly created object.

fill_log_normal(data, mean, stddev, stream=None)#

Fills in GPUArray data with log-normally distributed pseudorandom values with mean mean and standard deviation stddev.

CUDA 4.0 and above.

New in version 2012.2.

gen_log_normal(shape, dtype, mean, stddev, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with log-normally distributed pseudorandom values with mean mean and standard deviation stddev, and returns newly created object.

CUDA 4.0 and above.

New in version 2012.2.

fill_poisson(data, lambda_value=None, stream=None)#

Fills in GPUArray data with Poisson distributed pseudorandom values.

If lambda_value is not None, it is used as lambda, and data must be of type 32-bit unsigned int.

If lambda_value is None, the lambda value is read from each data array element (similarly to numpy.random.poisson), and the array is overwritten by the pseudorandom values. data must be of type 32-bit unsigned int, 32 or 64-bit float.

CUDA 5.0 and above.

New in version 2013.1.

gen_poisson(shape, dtype, lambda_value, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with Poisson distributed pseudorandom values with lambda lambda_value, and returns newly created object. dtype must be 32-bit unsigned int.

CUDA 5.0 and above.

New in version 2013.1.

call_skip_ahead(i, stream=None)#

Forces all generators to skip i values. Is equivalent to generating i values and discarding results, but is much faster.

call_skip_ahead_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many values to skip.

call_skip_ahead_sequence(i, stream=None)#

Forces all generators to skip i subsequences. Is equivalent to generating i * \(2^67\) values and discarding results, but is much faster.

call_skip_ahead_sequence_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many subsequences to skip.

class pycuda.curandom.MRG32k3aRandomNumberGenerator(seed_getter=None, offset=0)#
Parameters:
  • seed_getter – a function that, given an integer count, will yield an int32 GPUArray of seeds.

  • offset – Starting index into the XORWOW sequence, given seed.

Provides pseudorandom numbers. Generates sequences with period at least \(2^190\).

CUDA 4.1 and above.

New in version 2013.1.

fill_uniform(data, stream=None)#

Fills in GPUArray data with uniformly distributed pseudorandom values.

gen_uniform(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with uniformly distributed pseudorandom values, and returns newly created object.

fill_normal(data, stream=None)#

Fills in GPUArray data with normally distributed pseudorandom values.

gen_normal(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with normally distributed pseudorandom values, and returns newly created object.

fill_log_normal(data, mean, stddev, stream=None)#

Fills in GPUArray data with log-normally distributed pseudorandom values with mean mean and standard deviation stddev.

gen_log_normal(shape, dtype, mean, stddev, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with log-normally distributed pseudorandom values with mean mean and standard deviation stddev, and returns newly created object.

fill_poisson(data, lambda_value, stream=None)#

Fills in GPUArray data with Poisson distributed pseudorandom values.

If lambda_value is not None, it is used as lambda, and data must be of type 32-bit unsigned int.

If lambda_value is None, the lambda value is read from each data array element (similarly to numpy.random.poisson), and the array is overwritten by the pseudorandom values. data must be of type 32-bit unsigned int, 32 or 64-bit float.

CUDA 5.0 and above.

New in version 2013.1.

gen_poisson(shape, dtype, lambda_value, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with Poisson distributed pseudorandom values with lambda lambda_value, and returns newly created object. dtype must be 32-bit unsigned int.

CUDA 5.0 and above.

New in version 2013.1.

call_skip_ahead(i, stream=None)#

Forces all generators to skip i values. Is equivalent to generating i values and discarding results, but is much faster.

call_skip_ahead_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many values to skip.

call_skip_ahead_sequence(i, stream=None)#

Forces all generators to skip i subsequences. Is equivalent to generating i * \(2^67\) values and discarding results, but is much faster.

call_skip_ahead_sequence_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many subsequences to skip.

pycuda.curandom.generate_direction_vectors(count, direction=direction_vector_set.VECTOR_32)#

Return an GPUArray count filled with direction vectors used to initialize Sobol generators.

pycuda.curandom.generate_scramble_constants32(count)#

Return a GPUArray filled with count’ 32-bit unsigned integer numbers used to initialize :class:`ScrambledSobol32RandomNumberGenerator

pycuda.curandom.generate_scramble_constants64(count)#

Return a GPUArray filled with count’ 64-bit unsigned integer numbers used to initialize :class:`ScrambledSobol64RandomNumberGenerator

class pycuda.curandom.Sobol32RandomNumberGenerator(dir_vector=None, offset=0)#
Parameters:
  • dir_vector – a GPUArray of 32-element int32 vectors which are used to initialize quasirandom generator; it must contain one vector for each initialized generator

  • offset – Starting index into the Sobol32 sequence, given direction vector.

Provides quasirandom numbers. Generates sequences with period of \(2^32\).

CUDA 3.2 and above.

New in version 2011.1.

fill_uniform(data, stream=None)#

Fills in GPUArray data with uniformly distributed quasirandom values.

gen_uniform(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with uniformly distributed pseudorandom values, and returns newly created object.

fill_normal(data, stream=None)#

Fills in GPUArray data with normally distributed quasirandom values.

gen_normal(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with normally distributed pseudorandom values, and returns newly created object.

fill_log_normal(data, mean, stddev, stream=None)#

Fills in GPUArray data with log-normally distributed pseudorandom values with mean mean and standard deviation stddev.

CUDA 4.0 and above.

New in version 2012.2.

gen_log_normal(shape, dtype, mean, stddev, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with log-normally distributed pseudorandom values with mean mean and standard deviation stddev, and returns newly created object.

CUDA 4.0 and above.

New in version 2012.2.

fill_poisson(data, lambda_value, stream=None)#

Fills in GPUArray data with Poisson distributed pseudorandom values.

If lambda_value is not None, it is used as lambda, and data must be of type 32-bit unsigned int.

If lambda_value is None, the lambda value is read from each data array element (similarly to numpy.random.poisson), and the array is overwritten by the pseudorandom values. data must be of type 32-bit unsigned int, 32 or 64-bit float.

CUDA 5.0 and above.

New in version 2013.1.

gen_poisson(shape, dtype, lambda_value, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with Poisson distributed pseudorandom values with lambda lambda_value, and returns newly created object. dtype must be 32-bit unsigned int.

CUDA 5.0 and above.

New in version 2013.1.

call_skip_ahead(i, stream=None)#

Forces all generators to skip i values. Is equivalent to generating i values and discarding results, but is much faster.

call_skip_ahead_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many values to skip.

class pycuda.curandom.ScrambledSobol32RandomNumberGenerator(dir_vector=None, scramble_vector=None, offset=0)#
Parameters:
  • dir_vector – a GPUArray of 32-element uint32 vectors which are used to initialize quasirandom generator; it must contain one vector for each initialized generator

  • scramble_vector – a GPUArray of uint32 elements which are used to initialize quasirandom generator; it must contain one number for each initialized generator

  • offset – Starting index into the Sobol32 sequence, given direction vector.

Provides quasirandom numbers. Generates sequences with period of \(2^32\).

CUDA 4.0 and above.

New in version 2011.1.

fill_uniform(data, stream=None)#

Fills in GPUArray data with uniformly distributed quasirandom values.

gen_uniform(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with uniformly distributed pseudorandom values, and returns newly created object.

fill_normal(data, stream=None)#

Fills in GPUArray data with normally distributed quasirandom values.

gen_normal(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with normally distributed pseudorandom values, and returns newly created object.

fill_log_normal(data, mean, stddev, stream=None)#

Fills in GPUArray data with log-normally distributed pseudorandom values with mean mean and standard deviation stddev.

CUDA 4.0 and above.

New in version 2012.2.

gen_log_normal(shape, dtype, mean, stddev, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with log-normally distributed pseudorandom values with mean mean and standard deviation stddev, and returns newly created object.

CUDA 4.0 and above.

New in version 2012.2.

fill_poisson(data, lambda_value, stream=None)#

Fills in GPUArray data with Poisson distributed pseudorandom values.

If lambda_value is not None, it is used as lambda, and data must be of type 32-bit unsigned int.

If lambda_value is None, the lambda value is read from each data array element (similarly to numpy.random.poisson), and the array is overwritten by the pseudorandom values. data must be of type 32-bit unsigned int, 32 or 64-bit float.

CUDA 5.0 and above.

New in version 2013.1.

gen_poisson(shape, dtype, lambda_value, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with Poisson distributed pseudorandom values with lambda lambda_value, and returns newly created object. dtype must be 32-bit unsigned int.

CUDA 5.0 and above.

New in version 2013.1.

call_skip_ahead(i, stream=None)#

Forces all generators to skip i values. Is equivalent to generating i values and discarding results, but is much faster.

call_skip_ahead_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many values to skip.

class pycuda.curandom.Sobol64RandomNumberGenerator(dir_vector=None, offset=0)#
Parameters:
  • dir_vector – a GPUArray of 64-element uint64 vectors which are used to initialize quasirandom generator; it must contain one vector for each initialized generator

  • offset – Starting index into the Sobol64 sequence, given direction vector.

Provides quasirandom numbers. Generates sequences with period of \(2^64\).

CUDA 4.0 and above.

New in version 2011.1.

fill_uniform(data, stream=None)#

Fills in GPUArray data with uniformly distributed quasirandom values.

gen_uniform(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with uniformly distributed pseudorandom values, and returns newly created object.

fill_normal(data, stream=None)#

Fills in GPUArray data with normally distributed quasirandom values.

gen_normal(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with normally distributed pseudorandom values, and returns newly created object.

fill_log_normal(data, mean, stddev, stream=None)#

Fills in GPUArray data with log-normally distributed pseudorandom values with mean mean and standard deviation stddev.

CUDA 4.0 and above.

New in version 2012.2.

gen_log_normal(shape, dtype, mean, stddev, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with log-normally distributed pseudorandom values with mean mean and standard deviation stddev, and returns newly created object.

CUDA 4.0 and above.

New in version 2012.2.

fill_poisson(data, lambda_value, stream=None)#

Fills in GPUArray data with Poisson distributed pseudorandom values.

If lambda_value is not None, it is used as lambda, and data must be of type 32-bit unsigned int.

If lambda_value is None, the lambda value is read from each data array element (similarly to numpy.random.poisson), and the array is overwritten by the pseudorandom values. data must be of type 32-bit unsigned int, 32 or 64-bit float.

CUDA 5.0 and above.

New in version 2013.1.

gen_poisson(shape, dtype, lambda_value, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with Poisson distributed pseudorandom values with lambda lambda_value, and returns newly created object. dtype must be 32-bit unsigned int.

CUDA 5.0 and above.

New in version 2013.1.

call_skip_ahead(i, stream=None)#

Forces all generators to skip i values. Is equivalent to generating i values and discarding results, but is much faster.

call_skip_ahead_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many values to skip.

class pycuda.curandom.ScrambledSobol64RandomNumberGenerator(dir_vector=None, scramble_vector=None, offset=0)#
Parameters:
  • dir_vector – a GPUArray of 64-element uint64 vectors which are used to initialize quasirandom generator; it must contain one vector for each initialized generator

  • scramble_vector – a GPUArray of uint64 vectors which are used to initialize quasirandom generator; it must contain one vector for each initialized generator

  • offset – Starting index into the ScrambledSobol64 sequence, given direction vector.

Provides quasirandom numbers. Generates sequences with period of \(2^64\).

CUDA 4.0 and above.

New in version 2011.1.

fill_uniform(data, stream=None)#

Fills in GPUArray data with uniformly distributed quasirandom values.

gen_uniform(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with uniformly distributed pseudorandom values, and returns newly created object.

fill_normal(data, stream=None)#

Fills in GPUArray data with normally distributed quasirandom values.

gen_normal(shape, dtype, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with normally distributed pseudorandom values, and returns newly created object.

fill_log_normal(data, mean, stddev, stream=None)#

Fills in GPUArray data with log-normally distributed pseudorandom values with mean mean and standard deviation stddev.

CUDA 4.0 and above.

New in version 2012.2.

gen_log_normal(shape, dtype, mean, stddev, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with log-normally distributed pseudorandom values with mean mean and standard deviation stddev, and returns newly created object.

CUDA 4.0 and above.

New in version 2012.2.

fill_poisson(data, lambda_value, stream=None)#

Fills in GPUArray data with Poisson distributed pseudorandom values.

If lambda_value is not None, it is used as lambda, and data must be of type 32-bit unsigned int.

If lambda_value is None, the lambda value is read from each data array element (similarly to numpy.random.poisson), and the array is overwritten by the pseudorandom values. data must be of type 32-bit unsigned int, 32 or 64-bit float.

CUDA 5.0 and above.

New in version 2013.1.

gen_poisson(shape, dtype, lambda_value, stream=None)#

Creates object of GPUArray with given shape and dtype, fills it in with Poisson distributed pseudorandom values with lambda lambda_value, and returns newly created object. dtype must be 32-bit unsigned int.

CUDA 5.0 and above.

New in version 2013.1.

call_skip_ahead(i, stream=None)#

Forces all generators to skip i values. Is equivalent to generating i values and discarding results, but is much faster.

call_skip_ahead_array(i, stream=None)#

Accepts array i of integer values, telling each generator how many values to skip.

Single-pass Custom Expression Evaluation#

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.

class pycuda.elementwise.ElementwiseKernel(arguments, operation, name='kernel', keep=False, options=[], preamble='')#

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.

preamble specifies some source code that is included before the elementwise kernel specification. You may use this to include other files and/or define functions that are used by operation.

__call__(*args, range=None, slice=None)#

Invoke the generated scalar kernel. The arguments may either be scalars or GPUArray instances.

If range is given, it must be a slice object and specifies the range of indices i for which the operation is carried out.

If slice is given, it must be a slice object and specifies the range of indices i for which the operation is carried out, truncated to the container. Also, slice may contain negative indices to index relative to the end of the array.

If stream is given, it must be a pycuda.driver.Stream object, where the execution will be serialized.

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.)

Custom Reductions#

class pycuda.reduction.ReductionKernel(dtype_out, neutral, reduce_expr, map_expr=None, arguments=None, name='reduce_kernel', keep=False, options=[], preamble='', allocator=None)#

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. reduce_expr must be associative.

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.

__call__(*args, stream=None, out=None)#

Invoke the generated reduction kernel. The arguments may either be scalars or GPUArray instances. The reduction will be done on each entry of the first vector argument.

If stream is given, it must be a pycuda.driver.Stream object, where the execution will be serialized.

With out the resulting single-entry GPUArray can be specified. Because offsets are supported one can store results anywhere (e.g. out=a[3]).

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="x[i]*y[i]",
        arguments="float *x, float *y")

my_dot_prod = krnl(a, b).get()

Or by specifying the output:

from pycuda.curandom import rand as curand
a = curand((10, 200), dtype=np.float32)
red = ReductionKernel(np.float32, neutral=0,
                           reduce_expr="a+b",
                           arguments="float *in")
a_sum = gpuarray.empty(10, dtype=np.float32)
for i in range(10):
    red(a[i], out=a_sum[i])
assert(np.allclose(a_sum.get(), a.get().sum(axis=1)))

Parallel Scan / Prefix Sum#

class pycuda.scan.ExclusiveScanKernel(dtype, scan_expr, neutral, name_prefix='scan', options=[], preamble='')#

Generates a kernel that can compute a prefix sum using any associative operation given as scan_expr. scan_expr uses the formal values “a” and “b” to indicate two operands of an associative binary operation. neutral is the neutral element of scan_expr, obeying scan_expr(a, neutral) == a.

dtype specifies the type of the arrays being operated on. name_prefix is used for kernel names to ensure recognizability in profiles and logs. options is a list of compiler options to use when building. preamble specifies a string of code that is inserted before the actual kernels.

__call__(self, input_ary, output_ary=None, allocator=None, queue=None)#
class pycuda.scan.InclusiveScanKernel(dtype, scan_expr, neutral=None, name_prefix='scan', options=[], preamble='', devices=None)#

Works like ExclusiveScanKernel. Unlike the exclusive case, neutral is not required.

Here’s a usage example:

knl = InclusiveScanKernel(np.int32, "a+b")

n = 2**20-2**18+5
host_data = np.random.randint(0, 10, n).astype(np.int32)
dev_data = gpuarray.to_gpu(queue, host_data)

knl(dev_data)
assert (dev_data.get() == np.cumsum(host_data, axis=0)).all()

Custom data types in Reduction and Scan#

If you would like to use your own (struct/union/whatever) data types in scan and reduction, define those types in the preamble and let PyCUDA know about them using this function:

pycuda.tools.register_dtype(dtype, name)#

dtype is a numpy.dtype().

GPGPU Algorithms#

Bogdan Opanchuk’s reikna offers a variety of GPU-based algorithms (FFT, RNG, matrix multiplication) designed to work with pycuda.gpuarray.GPUArray objects.