dedalus.core.transforms

Spectral transform classes.

Module Contents

class BallRadialTransform(grid_shape, coeff_size, axis, ell_maps, regindex, regtotal, k, alpha, dtype=np.complex128, dealias_before_converting=None)

Abstract base class for all transforms.

backward(cdata, gdata, axis)
backward_reduced(cdata, gdata)
forward(gdata, cdata, axis)
forward_reduced(gdata, cdata)
N3c
N3g
alpha
dealias_before_converting = None
ell_maps
intertwiner
k
regindex
regtotal
class ComplexFFT(grid_size, coeff_size)

Abstract base class for complex-to-complex FFTs.

resize_coeffs(data_in, data_out, axis, rescale)

Resize and rescale coefficients in standard FFT format by intermediate padding/truncation.

class ComplexFourierMMT(grid_size, coeff_size)

Complex-to-complex Fourier MMT.

backward_matrix()

Build backward transform matrix.

forward_matrix()

Build forward transform matrix.

class ComplexFourierTransform(grid_size, coeff_size)

Abstract base class for complex-to-complex Fourier transforms.

Parameters:
  • grid_size (int) – Grid size (N) along transform dimension.

  • coeff_size (int) – Coefficient size (M) along transform dimension.

Notes

Let KN = (N - 1) // 2 be the maximum fully resolved (non-Nyquist) mode on the grid. Let KM = (M - 1) // 2 be the maximum retained mode in coeff space. Then K = min(KN, KM) is the maximum wavenumber used in the transforms. A unit-amplitude normalization is used.

Forward transform:
if abs(k) <= K:

F(k) = (1/N) sum_{x=0}^{N-1} f(x) exp(-2 pi i k x / N)

else:

F(k) = 0

Backward transform:

f(x) = sum_{k=-K}^{K} F(k) exp(2 pi i k x / N)

Coefficient ordering:

If M is odd, the ordering is [0, 1, 2, …, KM, KM+1, -KM, -KM+1, …, -1], where the Nyquist mode k = KM + 1 is zeroed in both directions. If M is even, the ordering is [0, 1, 2, …, KM, -KM, -KM+1, …, -1].

KM
KN
Kmax
M
N
property wavenumbers

One-dimensional global wavenumber array.

class CosineMMT(grid_size, coeff_size)

Cosine MMT.

backward_matrix()

Build backward transform matrix.

forward_matrix()

Build forward transform matrix.

class CosineTransform(grid_size, coeff_size)

Abstract base class for cosine transforms.

Parameters:
  • grid_size (int) – Grid size (N) along transform dimension.

  • coeff_size (int) – Coefficient size (M) along transform dimension.

Notes

Let KN = (N - 1) be the maximum (Nyquist) mode on the grid. Let KM = (M - 1) be the maximum retained mode in coeff space. Then K = min(KN, KM) is the maximum wavenumber used in the transforms. A unit-amplitude normalization is used.

Forward transform:
if k == 0:

a(k) = (1/N) sum_{x=0}^{N-1} f(x)

elif k <= K:

a(k) = (2/N) sum_{x=0}^{N-1} f(x) cos(pi k x / N)

else:

a(k) = 0

Backward transform:

f(x) = sum_{k=0}^{K} a(k) cos(pi k x / N)

Coefficient ordering:

The cosine coefficients are ordered simply as [a(0), a(1), a(2), …, a(KM)]

KM
KN
Kmax
M
N
property wavenumbers

One-dimensional global wavenumber array.

class DiskRadialTransform(grid_shape, basis_shape, axis, m_maps, s, k, alpha, dtype=np.complex128, dealias_before_converting=None)
backward_reduced(cdata, gdata)
forward_reduced(gdata, cdata)
Nmax
Nphi
alpha
dealias_before_converting = None
k
m_maps
s
class FFTPACKRealFFT(grid_size, coeff_size)

Real-to-real FFT using scipy.fftpack.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class FFTWBase(*args, rigor=None, **kw)

Abstract base class for FFTW transforms.

rigor = None
class FFTWComplexFFT(*args, rigor=None, **kw)

Complex-to-complex FFT using FFTW.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class FFTWDCT(*args, rigor=None, **kw)

Fast cosine transform using FFTW.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class FFTWFastChebyshevTransform(grid_size, coeff_size, a, b, a0, b0, **kw)

Fast ultraspherical transform using scipy.fft and spectral conversion.

class FFTWHalfComplexFFT(*args, rigor=None, **kw)

Real-to-real FFT using FFTW half-complex DFT.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

repack(cdata, temp, axis)

Repack into complex coefficients and rescale for unit-amplitude normalization.

unpack_rescale(temp, cdata, axis, rescale)

Unpack halfcomplex coefficients and rescale for unit-amplitude normalization.

class FFTWRealFFT(*args, rigor=None, **kw)

Real-to-real FFT using FFTW.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class FastChebyshevTransform(grid_size, coeff_size, a, b, a0, b0, **kw)

Abstract base class for fast Chebyshev transforms including ultraspherical conversion. Subclasses should inherit from this class, then a FastCosineTransform subclass.

KM_orig
Kmax_orig
M_orig
class FastCosineTransform(*args, **kw)

Abstract base class for fast cosine transforms.

resize_rescale_backward(data_in, data_out, axis, Kmax)

Resize by padding/trunction and rescale to unit amplitude.

resize_rescale_forward(data_in, data_out, axis, Kmax)

Resize by padding/trunction and rescale to unit amplitude.

backward_rescale_pos = 0.5
backward_rescale_zero = 1
forward_rescale_pos
forward_rescale_zero
class JacobiMMT(grid_size, coeff_size, a, b, a0, b0, dealias_before_converting=None)

Jacobi polynomial MMTs.

backward_matrix()

Build backward transform matrix.

forward_matrix()

Build forward transform matrix.

class JacobiTransform(grid_size, coeff_size, a, b, a0, b0, dealias_before_converting=None)

Abstract base class for Jacobi polynomial transforms.

Parameters:
  • grid_size (int) – Grid size (N) along transform dimension.

  • coeff_size (int) – Coefficient size (M) along transform dimension.

  • a (int) – Jacobi “a” parameter for polynomials.

  • b (int) – Jacobi “b” parameters for polynomials.

  • a0 (int) – Jacobi “a” parameter for the quadrature grid.

  • b0 (int) – Jacobi “b” parameter for the quadrature grid.

Notes

TODO: We need to define the normalization we use here.

M
N
a
a0
b
b0
dealias_before_converting = None
class NonSeparableTransform(grid_shape, coeff_size, axis, dtype)

Abstract base class for all transforms.

backward(cdata, gdata, axis)
forward(gdata, cdata, axis)
N2c
N2g
class PolynomialTransform(grid_size, coeff_size)

Abstract base class for all transforms.

backward(cdata, gdata, axis)
forward(gdata, cdata, axis)
static resize_reduced(data_in, data_out)

Resize data by padding/truncation.

class RealFFT(grid_size, coeff_size)

Abstract base class for real-to-real FFTs using real-to-complex algorithms.

repack_rescale(cdata, temp, axis, rescale)

Repack into complex coefficients and rescale for unit-amplitude normalization.

unpack_rescale(temp, cdata, axis, rescale)

Unpack complex coefficients and rescale for unit-amplitude normalization.

class RealFourierMMT(grid_size, coeff_size)

Real-to-real Fourier MMT.

backward_matrix()

Build backward transform matrix.

forward_matrix()

Build forward transform matrix.

class RealFourierTransform(grid_size, coeff_size)

Abstract base class for real-to-real Fourier transforms.

Parameters:
  • grid_size (int) – Grid size (N) along transform dimension.

  • coeff_size (int) – Coefficient size (M) along transform dimension.

Notes

Let KN = (N - 1) // 2 be the maximum fully resolved (non-Nyquist) mode on the grid. Let KM = (M - 1) // 2 be the maximum retained mode in coeff space. Then K = min(KN, KM) is the maximum wavenumber used in the transforms. A unit-amplitude normalization is used.

Forward transform:
if k == 0:

a(k) = (1/N) sum_{x=0}^{N-1} f(x) b(k) = 0

elif k <= K:

a(k) = (2/N) sum_{x=0}^{N-1} f(x) cos(-2 pi k x / N) b(k) = -(2/N) sum_{x=0}^{N-1} f(x) sin(-2 pi k x / N)

else:

a(k) = 0 b(k) = 0

Backward transform:

f(x) = sum_{k=0}^{K} a(k) cos(2 pi k x / N) - b(k) sin(2 pi k x / N)

Coefficient ordering:

The cosine and minus-sine coefficients are interleaved as [a(0), b(0), a(1), b(1), a(2), b(2), …, a(KM), b(KM)] where the k = 0 minus-sine mode is zeroed in both directions.

KM
KN
Kmax
M
N
property wavenumbers

One-dimensional global wavenumber array.

class SWSHColatitudeTransform(Ntheta, Lmax, m_maps, s)

Abstract base class for all transforms.

backward_reduced(cdata, gdata)
forward_reduced(gdata, cdata)
Lmax
Ntheta
m_maps
s
class ScipyComplexFFT(grid_size, coeff_size)

Complex-to-complex FFT using scipy.fft.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class ScipyDCT(*args, **kw)

Fast cosine transform using scipy.fft.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class ScipyFastChebyshevTransform(grid_size, coeff_size, a, b, a0, b0, **kw)

Fast ultraspherical transform using scipy.fft and spectral conversion.

class ScipyRealFFT(grid_size, coeff_size)

Real-to-real FFT using scipy.fft.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class SeparableMatrixTransform

Abstract base class for separable matrix-multiplication transforms.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

abstractmethod backward_matrix()

Build backward transform matrix.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

abstractmethod forward_matrix()

Build forward transform matrix.

class SeparableTransform

Abstract base class for transforms that only apply to one dimension, independent of all others.

abstractmethod backward(cdata, gdata, axis)

Apply backward transform along specified axis.

abstractmethod forward(gdata, cdata, axis)

Apply forward transform along specified axis.

class Transform

Abstract base class for all transforms.

backward_DFT(cdata, gdata, axis)
backward_disk(cdata, gdata, axis, k, s, local_m)

Apply bakward radial transform to data with fixed s and varying m.

forward_DFT(gdata, cdata, axis)
forward_disk(gdata, cdata, axis, k0, k, s, local_m)

Apply forward radial transform to data with fixed s and varying m.

reduce_array(data, axis)

Return reduced 3D view of array collapsed above and below specified axis.

reduced_view_3(data, axis)
reduced_view_4(data, axis)
reduced_view_5(data, axis)
register_transform(basis, name)

Decorator to add transform to basis class dictionary.

GET_DEALIAS_BEFORE_CONVERTING
GET_FFTW_RIGOR
logger