dedalus.core.transforms

Spectral transform classes.

Module Contents

logger
GET_FFTW_RIGOR
GET_DEALIAS_BEFORE_CONVERTING
register_transform(basis, name)

Decorator to add transform to basis class dictionary.

class Transform

Abstract base class for all transforms.

class SeparableTransform

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

abstract forward(gdata, cdata, axis)

Apply forward transform along specified axis.

abstract backward(cdata, gdata, axis)

Apply backward transform along specified axis.

class SeparableMatrixTransform

Abstract base class for separable matrix-multiplication transforms.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

abstract forward_matrix()

Build forward transform matrix.

abstract backward_matrix()

Build backward 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.

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

Jacobi polynomial MMTs.

forward_matrix()

Build forward transform matrix.

backward_matrix()

Build backward 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].

property wavenumbers

One-dimensional global wavenumber array.

class ComplexFourierMMT(grid_size, coeff_size)

Complex-to-complex Fourier MMT.

forward_matrix()

Build forward transform matrix.

backward_matrix()

Build backward transform matrix.

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 ScipyComplexFFT(grid_size, coeff_size)

Complex-to-complex FFT using scipy.fft.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

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

Abstract base class for FFTW transforms.

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

Complex-to-complex FFT using FFTW.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

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.

property wavenumbers

One-dimensional global wavenumber array.

class RealFourierMMT(grid_size, coeff_size)

Real-to-real Fourier MMT.

forward_matrix()

Build forward transform matrix.

backward_matrix()

Build backward transform matrix.

class FFTPACKRealFFT(grid_size, coeff_size)

Real-to-real FFT using scipy.fftpack.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

class RealFFT(grid_size, coeff_size)

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

unpack_rescale(temp, cdata, axis, rescale)

Unpack complex coefficients and rescale for unit-amplitude normalization.

repack_rescale(cdata, temp, axis, rescale)

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

class ScipyRealFFT(grid_size, coeff_size)

Real-to-real FFT using scipy.fft.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

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

Real-to-real FFT using FFTW.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

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

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

unpack_rescale(temp, cdata, axis, rescale)

Unpack halfcomplex coefficients and rescale for unit-amplitude normalization.

repack(cdata, temp, axis)

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

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

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

property wavenumbers

One-dimensional global wavenumber array.

class CosineMMT(grid_size, coeff_size)

Cosine MMT.

forward_matrix()

Build forward transform matrix.

backward_matrix()

Build backward transform matrix.

class FastCosineTransform(*args, **kw)

Abstract base class for fast cosine transforms.

resize_rescale_forward(data_in, data_out, axis, Kmax)

Resize by padding/trunction and rescale to unit amplitude.

resize_rescale_backward(data_in, data_out, axis, Kmax)

Resize by padding/trunction and rescale to unit amplitude.

class ScipyDCT(*args, **kw)

Fast cosine transform using scipy.fft.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward transform along specified axis.

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

Fast cosine transform using FFTW.

forward(gdata, cdata, axis)

Apply forward transform along specified axis.

backward(cdata, gdata, axis)

Apply backward 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.

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

Fast ultraspherical transform using scipy.fft and spectral conversion.

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

Fast ultraspherical transform using scipy.fft and spectral conversion.

reduced_view_3(data, axis)
reduced_view_4(data, axis)
reduced_view_5(data, axis)
class PolynomialTransform(grid_size, coeff_size)

Abstract base class for all transforms.

static resize_reduced(data_in, data_out)

Resize data by padding/truncation.

forward(gdata, cdata, axis)
backward(cdata, gdata, axis)
reduce_array(data, axis)

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

forward_DFT(gdata, cdata, axis)
backward_DFT(cdata, gdata, axis)
class NonSeparableTransform(grid_shape, coeff_size, axis, dtype)

Abstract base class for all transforms.

forward(gdata, cdata, axis)
backward(cdata, gdata, axis)
class SWSHColatitudeTransform(Ntheta, Lmax, m_maps, s)

Abstract base class for all transforms.

forward_reduced(gdata, cdata)
backward_reduced(cdata, gdata)
class DiskRadialTransform(grid_shape, basis_shape, axis, m_maps, s, k, alpha, dtype=np.complex128, dealias_before_converting=None)
forward_reduced(gdata, cdata)
backward_reduced(cdata, gdata)
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.

forward(gdata, cdata, axis)
backward(cdata, gdata, axis)
forward_reduced(gdata, cdata)
backward_reduced(cdata, gdata)
forward_disk(gdata, cdata, axis, k0, k, s, local_m)

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

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

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