dedalus.core.basis

Abstract and built-in classes for spectral bases.

Module Contents

logger
DEFAULT_LIBRARY
FFTW_RIGOR
DIRICHLET_PRECONDITIONING
class Basis

Base class for spectral bases.

These classes define methods for transforming, differentiating, and integrating corresponding series represented by their spectral coefficients.

Parameters
  • base_grid_size (int) – Number of grid points

  • interval (tuple of floats) – Spatial domain of basis

  • dealias (float, optional) – Fraction of modes to keep after dealiasing (default: 1.)

Variables
  • grid_dtype (dtype) – Grid data type

  • coeff_size (int) – Number of spectral coefficients

  • coeff_embed (int) – Padded number of spectral coefficients for transform

  • coeff_dtype (dtype) – Coefficient data type

property library
abstract set_dtype(grid_dtype)

Set transforms based on grid data type.

abstract forward(gdata, cdata, axis, meta, scale)

Grid-to-coefficient transform.

abstract backward(cdata, gdata, axis, meta, scale)

Coefficient-to-grid transform.

abstract differentiate(cdata, cderiv, axis)

Differentiate using coefficients.

abstract integrate(cdata, cint, axis)

Integrate over interval using coefficients.

abstract interpolate(cdata, cint, position, axis)

Interpolate in interval using coefficients.

grid_size(scale)

Compute scaled grid size.

grid_spacing(scale=1.0)

Compute grid spacings.

grid_array_object(domain, axis)

Grid array object.

grid_spacing_object(domain, axis)
check_arrays(cdata, gdata, axis, scale)

Verify provided arrays sizes and dtypes are correct. Build compliant arrays if not provided.

class TransverseBasis

Base class for bases supporting transverse differentiation.

class ImplicitBasis

Base class for bases supporting implicit methods.

These bases define the following matrices encoding the respective linear functions acting on a series represented by its spectral coefficients:

Linear operators (square matrices):

Pre : preconditioning (default: identity) Diff : differentiation Mult(p) : multiplication by p-th basis element

Linear functionals (vectors):

left_vector : left-endpoint evaluation right_vector : right-endpoint evaluation integ_vector : integration over interval interp_vector : interpolation in interval

abstract Multiply(p, ncc_basis_meta, arg_basis_meta)

p-element multiplication matrix.

abstract Precondition()

Preconditioning matrix.

DropTau(n_tau)

Matrix dropping tau+match rows.

DropNonfirst()

Matrix dropping non-first rows.

DropNonconstant()

Matrix dropping non-constant rows.

DropMatch()

Matrix dropping match rows.

PreconditionDropTau(n_tau)

Preconditioning with tau+match filtering.

PreconditionDropMatch()

Preconditioning with match filtering.

NCC(ncc_basis_meta, arg_basis_meta, coeffs, cutoff, max_terms)

Build NCC multiplication matrix.

class Chebyshev(name, base_grid_size, interval=(- 1, 1), dealias=1, tau_after_pre=True)

Chebyshev polynomial basis on the roots grid.

element_label = 'T'
separable = False
coupled = True
default_meta()
grid(scale=1.0)

Build Chebyshev roots grid.

grid_spacing(scale=1.0)

Build Chebyshev roots grid spacing.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

Integrate()

Build integration class.

Interpolate()

Buld interpolation class.

Differentiate()

Build differentiation class.

Precondition()

Preconditioning matrix.

T_n = (U_n - U_(n-2)) / 2 U_(-n) = -U_(n-2)

Dirichlet()

Dirichlet recombination matrix.

D[0] = T[0] D[1] = T[1] D[n] = T[n] - T[n-2]

<T[i]|D[j]> = <T[i]|T[j]> - <T[i]|T[j-2]>

= δ(i,j) - δ(i,j-2)

Multiply(p, ncc_basis_meta, arg_basis_meta)

p-element multiplication matrix.

class Legendre(name, base_grid_size, interval=(- 1, 1), dealias=1, tau_after_pre=True)

Legendre polynomial basis on the roots grid.

element_label = 'P'
separable = False
coupled = True
default_meta()
grid(scale=1.0)

Build Legendre polynomial roots grid.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

Integrate()

Build integration class.

Interpolate()

Buld interpolation class.

Differentiate()

Build differentiation class.

Precondition()
Preconditioning matrix.

2 * (2*n + 1) * P[n] = (n + 2) * J11[n] - n * J11[n-2]

Dirichlet()

Dirichlet recombination matrix.

D[0] = P[0] D[1] = P[1] D[n] = P[n] - P[n-2]

<P[i]|D[j]> = <P[i]|T[j]> - <P[i]|T[j-2]>

= δ(i,j) - δ(i,j-2)

Multiply(p, ncc_basis_meta, arg_basis_meta)

p-element multiplication matrix.

class Hermite(name, base_grid_size, center=0.0, stretch=1.0, dealias=1, tau_after_pre=False)

Hermite function/polynomial basis.

Interval:

The functions live on (-inf, inf) but are centered and scaled by the affine transformation from (-1, 1) to the specified interval.

Hermite functions (with envelope, default):

hf[n](xp) = exp(-xn^2/2) H[n](xn) / N[n] N[n]**2 = π^(1/2) 2^n n! integ hf[n] hf[m] dxn = δ[n,m]

Hermite polynomials (without envelope):

hp[n](xp) = H[n](xn) integ hp[n] hp[m] exp(-xn**2) dx = δ[n,m] N[n]**2

element_label = 'h'
separable = False
coupled = True
default_meta()
grid(scale=1.0)

Build Hermite polynomial roots grid.

grid_array_object(domain, axis)

Grid array object.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

Integrate()

Build integration class.

Interpolate()

Buld interpolation class.

Differentiate()

Build differentiation class.

Precondition()

Preconditioning matrix.

Dirichlet()

Dirichlet recombination matrix.

Multiply(p, ncc_basis_meta, arg_basis_meta)

Hermite multiplication matrix.

class Laguerre(name, base_grid_size, edge=0.0, stretch=1.0, dealias=1, tau_after_pre=False)

Laguerre function/polynomial basis.

Interval:

The functions live on (0, inf) but are centered and scaled by the affine transformation from (0, 1) to the specified interval.

Laguerre functions (with envelope, default):

gn[n](xp) = exp(-x/2) L[n](xn) integ gn[n] gn[m] dxn = δ[n,m]

Laguerre polynomials (without envelope):

gp[n](xp) = L[n](xn) integ gp[n] gp[m] exp(-xn) dx = δ[n,m]

element_label = 'g'
separable = False
coupled = True
default_meta()
grid(scale=1.0)

Build Laguerre polynomial roots grid.

grid_array_object(domain, axis)

Grid array object.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

Integrate()

Build integration class.

Interpolate()

Buld interpolation class.

Differentiate()

Build differentiation class.

Precondition()
Preconditioning matrix.

L[n;1] = sum_{i=0}^{n} L[n;0] L[n;0] = L[n;1] - L[n-1;1]

Dirichlet()

Dirichlet recombination matrix.

G[0] = g[0] G[n] = g[n] - g[n-1]

<g[i]|G[j]> = <g[i]|g[j]> - <g[i]|g[j-1]>

= δ(i,j) - δ(i,j-1)

Multiply(p, ncc_basis_meta, arg_basis_meta)

Laguerre multiplication matrix.

class Fourier(name, base_grid_size, interval=(0, 2 * pi), dealias=1)

Fourier complex exponential basis.

separable = True
coupled = False
element_label = 'k'
default_meta()
grid(scale=1.0)

Build evenly spaced Fourier grid.

grid_spacing(scale=1.0)

Build Fourier grid spacing.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

Integrate()

Build integration class.

Interpolate()

Build interpolation class.

Differentiate()

Build differentiation class.

HilbertTransform()

Build Hilbert transform class.

class SinCos(name, base_grid_size, interval=(0, pi), dealias=1)

Sin/Cos series basis.

element_label = 'k'
separable = True
coupled = False
default_meta()
grid(scale=1.0)

Evenly spaced interior grid: cos(Nx) = 0

grid_spacing(scale=1.0)

Build cosine grid spacing.

grid_array_object(domain, axis)

Grid array object.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

Integrate()

Build integration class.

Interpolate()

Build interpolation class.

Differentiate()

Build differentiation class.

HilbertTransform()

Build Hilbert transform class.

class Compound(name, subbases, dealias=1)

Compound basis joining adjascent subbases.

property library
separable = False
coupled = True
default_meta()
grid(scale=1.0)

Build compound grid.

grid_spacing(scale=1.0)

Build compound grid spacing.

grid_array_object(domain, axis)

Grid array object.

set_dtype(grid_dtype)

Determine coefficient properties from grid dtype.

coeff_start(index)
grid_start(index, scale)
sub_gdata(gdata, index, axis)

Retreive gdata corresponding to one subbasis.

sub_cdata(cdata, index, axis)

Retrieve cdata corresponding to one subbasis.

forward(gdata, cdata, axis, meta, scale)

Forward transforms.

backward(cdata, gdata, axis, meta, scale)

Backward transforms.

Integrate()

Build integration class.

Interpolate()

Buld interpolation class.

Differentiate()

Build differentiation class.

Precondition()

Preconditioning matrix.

Dirichlet()
Multiply(subindex, p, ncc_basis_meta, arg_basis_meta)

p-element multiplication matrix.

NCC(ncc_basis_meta, arg_basis_meta, coeffs, cutoff, max_terms)

Build NCC multiplication matrix.

DropTau(n_tau)

Matrix dropping tau+match rows.

DropNonfirst()

Matrix dropping non-first rows.

DropNonconstant()

Matrix dropping non-constant rows.

DropMatch()

Matrix dropping last row from each subbasis except the last.

PreconditionDropTau(n_tau)

Preconditioning and tau+match filtering.

PreconditionDropMatch()

Preconditioning and match filtering.

MatchRows()

Matrix matching subbases.