mpi4py_fft package¶
Subpackages¶
Submodules¶
mpi4py_fft.libfft module¶
-
class
mpi4py_fft.libfft.
FFT
(shape, axes=None, dtype=<class 'float'>, padding=False, backend='fftw', transforms=None, **kw)[source]¶ Bases:
mpi4py_fft.libfft.FFTBase
Class for serial FFT transforms
Parameters: shape (list or tuple of ints) – shape of input array planned for
axes (None, int or tuple of ints, optional) – axes to transform over. If None transform over all axes
dtype (np.dtype, optional) – Type of input array
padding (bool, number or list of numbers) – If False, then no padding. If number, then apply this number as padding factor for all axes. If list of numbers, then each number gives the padding for each axis. Must be same length as axes.
backend (str, optional) – Choose backend for serial transforms (
fftw
,pyfftw
,numpy
,scipy
,mkl_fft
). Default isfftw
transforms (None or dict, optional) – Dictionary of axes to serial transforms (forward and backward) along those axes. For example:
{(0, 1): (dctn, idctn), (2, 3): (dstn, idstn)}
If missing the default is to use rfftn/irfftn for real input arrays and fftn/ifftn for complex input arrays. Real-to-real transforms can be configured using this dictionary and real-to-real transforms from the
fftw.xfftn
module.kw (dict) – Parameters passed to serial transform object
-
forward
(input_array=None, output_array=None, **options)¶ Generic serial forward transform
Parameters: - input_array (array, optional)
- output_array (array, optional)
- options (dict) – parameters to serial transforms
Returns: output_array
Return type: array
-
backward
(input_array=None, output_array=None, **options)¶ Generic serial backward transform
Parameters: - input_array (array, optional)
- output_array (array, optional)
- options (dict) – parameters to serial transforms
Returns: output_array
Return type: array
-
class
mpi4py_fft.libfft.
FFTBase
(shape, axes=None, dtype=<class 'float'>, padding=False)[source]¶ Bases:
object
Base class for serial FFT transforms
Parameters: - shape (list or tuple of ints) – shape of input array planned for
- axes (None, int or tuple of ints, optional) – axes to transform over. If None transform over all axes
- dtype (np.dtype, optional) – Type of input array
- padding (bool, number or list of numbers) – If False, then no padding. If number, then apply this number as padding factor for all axes. If list of numbers, then each number gives the padding for each axis. Must be same length as axes.
mpi4py_fft.mpifft module¶
-
class
mpi4py_fft.mpifft.
PFFT
(comm, shape=None, axes=None, dtype=<class 'float'>, grid=None, padding=False, collapse=False, backend='fftw', transforms=None, darray=None, **kw)[source]¶ Bases:
object
Base class for parallel FFT transforms
Parameters: comm (MPI communicator)
shape (sequence of ints, optional) – shape of input array planned for
axes (None, int, sequence of ints or sequence of sequence of ints, optional) – axes to transform over.
- None -> All axes are transformed
- int -> Just one axis to transform over
- sequence of ints -> e.g., (0, 1, 2) or (0, 2, 1)
- sequence of sequence of ints -> e.g., ((0,), (1,)) or ((0,), (1, 2)) For seq. of seq. of ints all but the last transformed sequence may be longer than 1. This corresponds to collapsing axes, where serial FFTs are performed for all collapsed axes in one single call
dtype (np.dtype, optional) – Type of input array
grid (sequence of ints, optional) – Define processor grid sizes. Non positive values act as wildcards to allow MPI compute optimal decompositions. The sequence is padded with ones to match the global transform dimension. Use
(-1,)
to get a slab decomposition on the first axis. Use(1, -1)
to get a slab decomposition on the second axis. Use(P, Q)
or(P, Q, 1)
to get a 3D transform with 2D-pencil decomposition on a PxQ processor grid with the last axis non distributed. Use(P, 1, Q)
to get a 3D transform with 2D-pencil decomposition on a PxQ processor grid with the second to last axis non distributed.padding (bool, number or sequence of numbers, optional) – If False, then no padding. If number, then apply this number as padding factor for all axes. If sequence of numbers, then each number gives the padding for each axis. Must be same length as axes.
collapse (bool, optional) – If True try to collapse several serial transforms into one
backend (str, optional) – Choose backend for serial transforms (
fftw
,pyfftw
,numpy
,scipy
,mkl_fft
). Default isfftw
transforms (None or dict, optional) – Dictionary of axes to serial transforms (forward and backward) along those axes. For example:
{(0, 1): (dctn, idctn), (2, 3): (dstn, idstn)}
If missing the default is to use rfftn/irfftn for real input arrays and fftn/ifftn for complex input arrays. Real-to-real transforms can be configured using this dictionary and real-to-real transforms from the
fftw.xfftn
module. See Examples.
Other Parameters: - darray (DistArray object, optional) – Create PFFT using information contained in
darray
, neglecting most optional Parameters above - slab (bool or int, optional) – DEPRECATED. If True then distribute only one axis of the global array.
-
forward
(input_array=None, output_array=None, **kw)¶ Parallel forward transform. The method is an instance of the
Transform
class. SeeTransform.__call__()
Parameters: - input_array (array, optional)
- output_array (array, optional)
- kw (dict) – parameters to serial transforms
Returns: output_array
Return type: array
-
backward
(input_array=None, output_array=None, **kw)¶ Parallel backward transform. The method is an instance of the
Transform
class. SeeTransform.__call__()
Parameters: - input_array (array, optional)
- output_array (array, optional)
- kw (dict) – parameters to serial transforms
Returns: output_array
Return type: array
Examples
>>> import numpy as np >>> from mpi4py import MPI >>> from mpi4py_fft import PFFT, newDistArray >>> N = np.array([12, 14, 15], dtype=int) >>> fft = PFFT(MPI.COMM_WORLD, N, axes=(0, 1, 2)) >>> u = newDistArray(fft, False) >>> u[:] = np.random.random(u.shape).astype(u.dtype) >>> u_hat = fft.forward(u) >>> uj = np.zeros_like(u) >>> uj = fft.backward(u_hat, uj) >>> assert np.allclose(uj, u)
Now configure with real-to-real discrete cosine transform type 3
>>> from mpi4py_fft.fftw import rfftn, irfftn, dctn, idctn >>> import functools >>> dct = functools.partial(dctn, type=3) >>> idct = functools.partial(idctn, type=3) >>> transforms = {(1, 2): (dct, idct)} >>> r2c = PFFT(MPI.COMM_WORLD, N, axes=((0,), (1, 2)), transforms=transforms) >>> u = newDistArray(r2c, False) >>> u[:] = np.random.random(u.shape).astype(u.dtype) >>> u_hat = r2c.forward(u) >>> uj = np.zeros_like(u) >>> uj = r2c.backward(u_hat, uj) >>> assert np.allclose(uj, u)
-
dimensions
¶ The number of dimensions for transformed arrays
-
dtype
(forward_output=False)[source]¶ The type of transformed arrays
Parameters: forward_output (bool, optional) – If True then return dtype of an array that is the result of a forward transform. Otherwise, return the dtype of an array that is input to a forward transform.
-
global_shape
(forward_output=False)[source]¶ Return global shape of associated tensors
Parameters: forward_output (bool, optional) – If True then return global shape of spectral space, i.e., the input to a backward transfer. If False then return shape of physical space, i.e., the input to a forward transfer.
-
class
mpi4py_fft.mpifft.
Transform
(xfftn, transfer, pencil)[source]¶ Bases:
object
Class for performing any parallel transform, forward or backward
Parameters: - xfftn (list of serial transform objects)
- transfer (list of global redistribution objects)
- pencil (list of two pencil objects) – The two pencils represent the input and final output configuration of the distributed global arrays
-
__call__
(input_array=None, output_array=None, **kw)[source]¶ Compute transform
Parameters: - input_array (array, optional)
- output_array (array, optional)
- kw (dict) – parameters to serial transforms Note in particular that the keyword ‘normalize’=True/False can be used to turn normalization on or off. Default is to enable normalization for forward transforms and disable it for backward.
Note
If input_array/output_array are not given, then use predefined arrays as planned with serial transform object _xfftn.
-
input_array
¶ Return input array of Transform
-
input_pencil
¶ Return input pencil of Transform
-
output_array
¶ Return output array of Transform
-
output_pencil
¶ Return output pencil of Transform
mpi4py_fft.pencil module¶
-
class
mpi4py_fft.pencil.
Pencil
(subcomm, shape, axis=-1)[source]¶ Bases:
object
Class to represent a distributed array (pencil)
Parameters: - subcomm (MPI communicator)
- shape (sequence of ints) – Shape of global array
- axis (int, optional) – Pencil is aligned in this direction
Examples
Create two pencils for a 4-dimensional array of shape (8, 8, 8, 8) using 4 processors in total. The input pencil will be distributed in the first two axes, whereas the output pencil will be distributed in axes 1 and 2. Note that the Subcomm instance below may distribute any axis where an entry 0 is found, whereas an entry of 1 means that this axis should not be distributed.
>>> import subprocess >>> fx = open('pencil_script.py', 'w') >>> h = fx.write(''' ... import numpy as np ... from mpi4py import MPI ... from mpi4py_fft.pencil import Subcomm, Pencil ... comm = MPI.COMM_WORLD ... N = (8, 8, 8, 8) ... subcomms = Subcomm(comm, [0, 0, 1, 0]) ... axis = 2 ... p0 = Pencil(subcomms, N, axis) ... p1 = p0.pencil(0) ... shape0 = comm.gather(p0.subshape) ... shape1 = comm.gather(p1.subshape) ... if comm.Get_rank() == 0: ... print('Subshapes all 4 processors pencil p0:') ... print(np.array(shape0)) ... print('Subshapes all 4 processors pencil p1:') ... print(np.array(shape1))''') >>> fx.close() >>> print(subprocess.getoutput('mpirun -np 4 python pencil_script.py')) Subshapes all 4 processors pencil p0: [[4 4 8 8] [4 4 8 8] [4 4 8 8] [4 4 8 8]] Subshapes all 4 processors pencil p1: [[8 4 4 8] [8 4 4 8] [8 4 4 8] [8 4 4 8]]
Two index sets of the global data of shape (8, 8, 8, 8) are distributed. This means that the current distribution is using two groups of processors, with 2 processors in each group (4 in total). One group shares axis 0 and the other axis 1 on the input arrays. On the output, one group shares axis 1, whereas the other shares axis 2. Note that the call
p1 = p0.pencil(0)
creates a new pencil (p1) that is non-distributed in axes 0. It is, in other words, aligned in axis 0. Hence the first 8 in the lists with [8 4 4 8] above. The alignment is configurable, andp1 = p0.pencil(1)
would lead to an output pencil aligned in axis 1.-
pencil
(axis)[source]¶ Return a Pencil aligned with axis
Parameters: axis (int) – The axis along which the pencil is aligned
-
transfer
(pencil, dtype)[source]¶ Return an appropriate instance of the
Transfer
classThe returned
Transfer
class is used for global redistribution from this pencil’s instance to the pencil instance provided.Parameters: - pencil (
Pencil
) – The receiving pencil of a forward transform - dtype (dtype) – The type of the sending pencil
- pencil (
-
class
mpi4py_fft.pencil.
Subcomm
[source]¶ Bases:
tuple
Class returning a tuple of subcommunicators of any dimensionality
Parameters: - comm (A communicator or group of communicators)
- dims (None, int or sequence of ints) – dims = [0, 0, 1] will give communicators distributed in the two first indices, whereas the third will not be distributed
Examples
>>> import subprocess >>> fx = open('subcomm_script.py', 'w') >>> h = fx.write(''' ... from mpi4py import MPI ... comm = MPI.COMM_WORLD ... from mpi4py_fft.pencil import Subcomm ... subcomms = Subcomm(comm, [0, 0, 1]) ... if comm.Get_rank() == 0: ... for subcomm in subcomms: ... print(subcomm.Get_size())''') >>> fx.close() >>> print(subprocess.getoutput('mpirun -np 4 python subcomm_script.py')) 2 2 1 >>> print(subprocess.getoutput('mpirun -np 6 python subcomm_script.py')) 3 2 1
-
class
mpi4py_fft.pencil.
Transfer
(comm, shape, dtype, subshapeA, axisA, subshapeB, axisB)[source]¶ Bases:
object
Class for performing global redistributions
Parameters: - comm (MPI communicator)
- shape (sequence of ints) – shape of input array planned for
- dtype (np.dtype, optional) – Type of input array
- subshapeA (sequence of ints) – Shape of input pencil
- axisA (int) – Input array aligned in this direction
- subshapeB (sequence of ints) – Shape of output pencil
- axisB (int) – Output array aligned in this direction
Examples
Create two pencils for a 4-dimensional array of shape (8, 8, 8, 8) using 4 processors in total. The input pencil will be distributed in the first two axes, whereas the output pencil will be distributed in axes 1 and 2. Create a random array of shape according to the input pencil and transfer its values to an array of the output shape.
>>> import subprocess >>> fx = open('transfer_script.py', 'w') >>> h = fx.write(''' ... import numpy as np ... from mpi4py import MPI ... from mpi4py_fft.pencil import Subcomm, Pencil ... comm = MPI.COMM_WORLD ... N = (8, 8, 8, 8) ... subcomms = Subcomm(comm, [0, 0, 1, 0]) ... axis = 2 ... p0 = Pencil(subcomms, N, axis) ... p1 = p0.pencil(0) ... transfer = p0.transfer(p1, np.float) ... a0 = np.zeros(p0.subshape, dtype=np.float) ... a1 = np.zeros(p1.subshape) ... a0[:] = np.random.random(a0.shape) ... transfer.forward(a0, a1) ... s0 = comm.reduce(np.sum(a0**2)) ... s1 = comm.reduce(np.sum(a1**2)) ... if comm.Get_rank() == 0: ... assert np.allclose(s0, s1)''') >>> fx.close() >>> h=subprocess.getoutput('mpirun -np 4 python transfer_script.py')
mpi4py_fft.distarray module¶
-
class
mpi4py_fft.distarray.
DistArray
[source]¶ Bases:
numpy.ndarray
Distributed Numpy array
This Numpy array is part of a larger global array. Information about the distribution is contained in the attributes.
Parameters: - global_shape (sequence of ints) – Shape of non-distributed global array
- subcomm (None,
Subcomm
object or sequence of ints, optional) – Describes how to distribute the array - val (Number or None, optional) – Initialize array with this number if buffer is not given
- dtype (np.dtype, optional) – Type of array
- buffer (Numpy array, optional) – Array of correct shape. The buffer owns the memory that is used for this array.
- alignment (None or int, optional) – Make sure array is aligned in this direction. Note that alignment does not take rank into consideration.
- rank (int, optional) – Rank of tensor (number of free indices, a scalar is zero, vector one, matrix two)
For more information, see numpy.ndarray
Note
Tensors of rank higher than 0 are not distributed in the first
rank
indices. For example,>>> from mpi4py_fft import DistArray >>> a = DistArray((3, 8, 8, 8), rank=1) >>> print(a.pencil.shape) (8, 8, 8)
The array
a
cannot be distributed in the first axis of length 3 since rank is 1 and this first index represent the vector component. Thepencil
attribute ofa
thus only considers the last three axes.Also note that the
alignment
keyword does not take rank into consideration. Setting alignment=2 for the array above means that the last axis will be aligned, also when rank>0.-
alignment
¶ Return alignment of local
self
arrayNote
For tensors of rank > 0 the array is actually aligned along
alignment+rank
-
commsizes
¶ Return number of processors along each axis of
self
-
dimensions
¶ Return dimensions of array not including rank
-
get
(gslice)[source]¶ Return global slice of
self
Parameters: gslice (sequence of slice(None) and ints) – The slice of the global array. Returns: The slice of the global array is returned on rank 0, whereas the remaining ranks return None Return type: Numpy array Example
>>> import subprocess >>> fx = open('gs_script.py', 'w') >>> h = fx.write(''' ... from mpi4py import MPI ... from mpi4py_fft.distarray import DistArray ... comm = MPI.COMM_WORLD ... N = (6, 6, 6) ... z = DistArray(N, dtype=float, alignment=0) ... z[:] = comm.Get_rank() ... g = z.get((0, slice(None), 0)) ... if comm.Get_rank() == 0: ... print(g)''') >>> fx.close() >>> print(subprocess.getoutput('mpirun -np 4 python gs_script.py')) [0. 0. 0. 2. 2. 2.]
-
get_pencil_and_transfer
(axis)[source]¶ Return pencil and transfer objects for alignment along
axis
Parameters: axis (int) – The new axis to align data with Returns: 2-tuple where first item is a Pencil
aligned inaxis
. Second item is aTransfer
object for executing the redistribution of dataReturn type: 2-tuple
-
global_shape
¶ Return global shape of
self
-
local_slice
()[source]¶ Return local view into global
self
arrayReturns: Each item of the returned list is the slice along that axis, describing the view of the self
array into the global array.Return type: List of slices Example
Print local_slice of a global array of shape (16, 14, 12) using 4 processors.
>>> import subprocess >>> fx = open('ls_script.py', 'w') >>> h = fx.write(''' ... from mpi4py import MPI ... from mpi4py_fft.distarray import DistArray ... comm = MPI.COMM_WORLD ... N = (16, 14, 12) ... z = DistArray(N, dtype=float, alignment=0) ... ls = comm.gather(z.local_slice()) ... if comm.Get_rank() == 0: ... for l in ls: ... print(l)''') >>> fx.close() >>> print(subprocess.getoutput('mpirun -np 4 python ls_script.py')) (slice(0, 16, None), slice(0, 7, None), slice(0, 6, None)) (slice(0, 16, None), slice(0, 7, None), slice(6, 12, None)) (slice(0, 16, None), slice(7, 14, None), slice(0, 6, None)) (slice(0, 16, None), slice(7, 14, None), slice(6, 12, None))
-
pencil
¶ Return pencil describing distribution of
self
-
rank
¶ Return tensor rank of
self
-
read
(filename, name='darray', step=0)[source]¶ Read data
name
at indexstep``from file ``filename
intoself
Note
Only whole arrays can be read from file, not slices.
Parameters: - filename (str or instance of
FileBase
) – The name of the file (or the file itself) holding the data that is loaded intoself
. - name (str, optional) – Internal name in file of snapshot to be read.
- step (int, optional) – Index of field to be read. Default is 0.
Example
>>> from mpi4py_fft import DistArray >>> u = DistArray((8, 8), val=1) >>> u.write('h5file.h5', 'u', 0) >>> v = DistArray((8, 8)) >>> v.read('h5file.h5', 'u', 0) >>> assert np.allclose(u, v)
- filename (str or instance of
-
redistribute
(axis=None, out=None)[source]¶ Global redistribution of local
self
arrayParameters: - axis (int, optional) – Align local
self
array along this axis - out (
DistArray
, optional) – Copy data to this array of possibly different alignment
Returns: DistArray – The
self
array globally redistributed. If keywordout
is None then a new DistArray (aligned alongaxis
) is created and returned. Otherwise the provided out array is returned.Return type: out
- axis (int, optional) – Align local
-
subcomm
¶ Return tuple of subcommunicators for all axes of
self
-
substart
¶ Return starting indices of local
self
array
-
v
¶ Return local
self
array as anndarray
object
-
write
(filename, name='darray', step=0, global_slice=None, domain=None, as_scalar=False)[source]¶ Write snapshot
step
ofself
to filefilename
Parameters: filename (str or instance of
FileBase
) – The name of the file (or the file itself) that is used to store the requested data inself
name (str, optional) – Name used for storing snapshot in file.
step (int, optional) – Index used for snapshot in file.
global_slice (sequence of slices or integers, optional) – Store only this global slice of
self
domain (sequence, optional) – An optional spatial mesh or domain to go with the data. Sequence of either
- 2-tuples, where each 2-tuple contains the (origin, length) of each dimension, e.g., (0, 2*pi).
- Arrays of coordinates, e.g., np.linspace(0, 2*pi, N). One array per dimension
as_scalar (boolean, optional) – Whether to store rank > 0 arrays as scalars. Default is False.
Example
>>> from mpi4py_fft import DistArray >>> u = DistArray((8, 8), val=1) >>> u.write('h5file.h5', 'u', 0) >>> u.write('h5file.h5', 'u', (slice(None), 4))
-
mpi4py_fft.distarray.
newDistArray
(pfft, forward_output=True, val=0, rank=0, view=False)[source]¶ Return a new
DistArray
object for providedPFFT
objectParameters: - pfft (
PFFT
object) - forward_output (boolean, optional) – If False then create DistArray of shape/type for input to forward transform, otherwise create DistArray of shape/type for output from forward transform.
- val (int or float, optional) – Value used to initialize array.
- rank (int, optional) – Scalar has rank 0, vector 1 and matrix 2.
- view (bool, optional) – If True return view of the underlying Numpy array, i.e., return cls.view(np.ndarray). Note that the DistArray still will be accessible through the base attribute of the view.
Returns: A new
DistArray
object. Return thendarray
view if keywordview
is True.Return type: Examples
>>> from mpi4py import MPI >>> from mpi4py_fft import PFFT, newDistArray >>> FFT = PFFT(MPI.COMM_WORLD, [64, 64, 64]) >>> u = newDistArray(FFT, False, rank=1) >>> u_hat = newDistArray(FFT, True, rank=1)
- pfft (
Module contents¶
This is the mpi4py-fft package
What is mpi4py-fft?¶
The Python package mpi4py-fft is a tool primarily for working with Fast Fourier Transforms (FFTs) of (large) multidimensional arrays. There is really no limit as to how large the arrays can be, just as long as there is sufficient computing powers available. Also, there are no limits as to how transforms can be configured. Just about any combination of transforms from the FFTW library is supported. Furthermore, mpi4py-fft can also be used simply to perform global redistributions (distribute and communicate) of large arrays with MPI, without any transforms at all.
For more information, see documentation.