API Reference

Main interface

class dtcwt.Pyramid(lowpass, highpasses, scales=None)

A representation of a transform domain signal.

Backends are free to implement any class which respects this interface for storing transform-domain signals. The inverse transform may accept a backend-specific version of this class but should always accept any class which corresponds to this interface.

lowpass

A NumPy-compatible array containing the coarsest scale lowpass signal.

highpasses

A tuple where each element is the complex subband coefficients for corresponding scales finest to coarsest.

scales

(optional) A tuple where each element is a NumPy-compatible array containing the lowpass signal for corresponding scales finest to coarsest. This is not required for the inverse and may be None.

class dtcwt.Transform1d(biort='near_sym_a', qshift='qshift_a')

An implementation of the 1D DT-CWT in NumPy.

Parameters
forward(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT decompostion on a 1D column vector X (or on the columns of a matrix X).

Parameters
  • X – 1D real array or 2D real array whose columns are to be transformed

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.Pyramid-like object representing the transform result.

If biort or qshift are strings, they are used as an argument to the biort() or qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

inverse(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 1D reconstruction.

Parameters
  • pyramid – A dtcwt.Pyramid-like object containing the transformed signal.

  • gain_mask – Gain to be applied to each subband.

Returns

Reconstructed real array.

The l-th element of gain_mask is gain for wavelet subband at level l. If gain_mask[l] == 0, no computation is performed for band l. Default gain_mask is all ones. Note that l is 0-indexed.

class dtcwt.Transform2d(biort='near_sym_a', qshift='qshift_a')

An implementation of the 2D DT-CWT via NumPy. biort and qshift are the wavelets which parameterise the transform.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

forward(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT-2D decompostion on a 2D matrix X.

Parameters
  • X – 2D real array

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.Pyramid compatible object representing the transform-domain signal

inverse(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 2D reconstruction.

Parameters
  • pyramid – A dtcwt.Pyramid-like class holding the transform domain representation to invert.

  • gain_mask – Gain to be applied to each subband.

Returns

A numpy-array compatible instance with the reconstruction.

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

class dtcwt.Transform3d(biort='near_sym_a', qshift='qshift_a', ext_mode=4)

An implementation of the 3D DT-CWT via NumPy. biort and qshift are the wavelets which parameterise the transform. Valid values are documented in dtcwt.coeffs.biort() and dtcwt.coeffs.qshift().

forward(X, nlevels=3, include_scale=False, discard_level_1=False)

Perform a n-level DTCWT-3D decompostion on a 3D matrix X.

Parameters
  • X – 3D real array-like object

  • nlevels – Number of levels of wavelet decomposition

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • discard_level_1 – True if level 1 high-pass bands are to be discarded.

Returns

a dtcwt.Pyramid instance

Each element of the Pyramid highpasses tuple is a 4D complex array with the 4th dimension having size 28. The 3D slice [l][:,:,:,d] corresponds to the complex higpass coefficients for direction d at level l where d and l are both 0-indexed.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

There are two values for ext_mode, either 4 or 8. If ext_mode = 4, check whether 1st level is divisible by 2 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coefs can be divided by 4. If any dimension size is not a multiple of 4, append extra coefs by repeating the edges. If ext_mode = 8, check whether 1st level is divisible by 4 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coeffs can be divided by 8. If any dimension size is not a multiple of 8, append extra coeffs by repeating the edges twice.

If discard_level_1 is True the highpass coefficients at level 1 will not be discarded. (And, in fact, will never be calculated.) This turns the transform from being 8:1 redundant to being 1:1 redundant at the cost of no-longer allowing perfect reconstruction. If this option is selected then the first element of the highpasses tuple will be None. Note that dtcwt.Transform3d.inverse() will accept the first element being None and will treat it as being zero.

inverse(pyramid)

Perform an n-level dual-tree complex wavelet (DTCWT) 3D reconstruction.

Parameters
  • pyramid – The dtcwt.Pyramid-like instance representing the transformed signal.

  • biort – Level 1 wavelets to use. See biort().

  • qshift – Level >= 2 wavelets to use. See qshift().

  • ext_mode – Extension mode. See below.

Returns

Reconstructed real image matrix.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

There are two values for ext_mode, either 4 or 8. If ext_mode = 4, check whether 1st level is divisible by 2 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coefs can be divided by 4. If any dimension size is not a multiple of 4, append extra coefs by repeating the edges. If ext_mode = 8, check whether 1st level is divisible by 4 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coeffs can be divided by 8. If any dimension size is not a multiple of 8, append extra coeffs by repeating the edges twice.

dtcwt.backend_name = 'numpy'

A string providing a short human-readable name for the DTCWT backend currently being used. This corresponds to the name parameter passed to dtcwt.push_backend(). The default backend is numpy but can be overridden by setting the DTCWT_BACKEND environment variable to a valid backend name.

dtcwt.pop_backend()

Restore the backend after a call to push_backend(). Calls to pop_backend() and pop_backend() may be nested. This function will undo the most recent call to push_backend().

Raises

IndexError – if one attempts to pop more backends than one has pushed.

dtcwt.preserve_backend_stack()

Return a generator object which can be used to preserve the backend stack even when an exception has been raise. For example:

# current backend is NumPy
assert dtcwt.backend_name == 'numpy'

with dtcwt.preserve_backend_stack():
    dtcwt.push_backend('opencl')
    # ... things which may raise an exception

# current backend is NumPy even if an exception was thrown
assert dtcwt.backend_name == 'numpy'
dtcwt.push_backend(name)

Switch backend implementation to name. Push the previous backend onto the backend stack. The previous backend may be restored via dtcwt.pop_backend().

Parameters

name – string identifying which backend to switch to

Raises

ValueError – if name does not correspond to a known backend

name may take one of the following values:

  • numpy: the default NumPy backend. See dtcwt.numpy.

  • opencl: a backend which uses OpenCL where available. See dtcwt.opencl.

Functions to load standard wavelet coefficients.

dtcwt.coeffs.biort(name)

Load level 1 wavelet by name.

Parameters

name – a string specifying the wavelet family name

Returns

a tuple of vectors giving filter coefficients

Name

Wavelet

antonini

Antonini 9,7 tap filters.

legall

LeGall 5,3 tap filters.

near_sym_a

Near-Symmetric 5,7 tap filters.

near_sym_b

Near-Symmetric 13,19 tap filters.

near_sym_b_bp

Near-Symmetric 13,19 tap filters + BP filter

Return a tuple whose elements are a vector specifying the h0o, g0o, h1o and g1o coefficients.

See Rotational symmetry modified wavelet transform for an explanation of the near_sym_b_bp wavelet filters.

Raises
  • IOError – if name does not correspond to a set of wavelets known to the library.

  • ValueError – if name specifies a dtcwt.coeffs.qshift() wavelet.

dtcwt.coeffs.qshift(name)

Load level >=2 wavelet by name,

Parameters

name – a string specifying the wavelet family name

Returns

a tuple of vectors giving filter coefficients

Name

Wavelet

qshift_06

Quarter Sample Shift Orthogonal (Q-Shift) 10,10 tap filters, (only 6,6 non-zero taps).

qshift_a

Q-shift 10,10 tap filters, (with 10,10 non-zero taps, unlike qshift_06).

qshift_b

Q-Shift 14,14 tap filters.

qshift_c

Q-Shift 16,16 tap filters.

qshift_d

Q-Shift 18,18 tap filters.

qshift_b_bp

Q-Shift 18,18 tap filters + BP

Return a tuple whose elements are a vector specifying the h0a, h0b, g0a, g0b, h1a, h1b, g1a and g1b coefficients.

See Rotational symmetry modified wavelet transform for an explanation of the qshift_b_bp wavelet filters.

Raises
  • IOError – if name does not correspond to a set of wavelets known to the library.

  • ValueError – if name specifies a dtcwt.coeffs.biort() wavelet.

Keypoint analysis

dtcwt.keypoint.find_keypoints(highpass_highpasses, method=None, alpha=1.0, beta=0.4, kappa=0.16666666666666666, threshold=None, max_points=None, upsample_keypoint_energy=None, upsample_highpasses=None, refine_positions=True, skip_levels=1)
Parameters
  • highpass_highpasses – (NxMx6) matrix of highpass subband images

  • method(optional) string specifying which keypoint energy method to use

  • alpha(optional) scale parameter for 'fauqueur' method

  • beta(optional) shape parameter for 'fauqueur' method

  • kappa(optiona) suppression parameter for 'kingsbury' method

  • threshold(optional) minimum keypoint energy of returned keypoints

  • max_points(optional) maximum number of keypoints to return

  • upsample_keypoint_energy – is non-None, a string specifying a method used to upscale the keypoint energy map before finding keypoints

  • upsample_subands – is non-None, a string specifying a method used to upscale the subband image before finding keypoints

  • refine_positions(optional) should the keypoint positions be refined to sub-pixel accuracy

  • skip_levels(optional) number of levels of the transform to ignore before looking for keypoints

Returns

(Px4) array of P keypoints in image co-ordinates

Warning

The interface and behaviour of this function is the subject of an open research project. It is provided in this release as a preview of forthcoming functionality but it is subject to change between releases.

The rows of the returned keypoint array give the x co-ordinate, y co-ordinate, scale and keypoint energy. The rows are sorted in order of decreasing keypoint energy.

If refine_positions is True then the positions (and energy) of the keypoints will be refined to sub-pixel accuracy by fitting a quadratic patch. If refine_positions is False then the keypoint locations will be those corresponding directly to pixel-wise maxima of the subband images.

The max_points and threshold parameters are cumulative: if both are specified then the max_points greatest energy keypoints with energy greater than threshold will be returned.

Usually the keypoint energies returned from the finest scale level are dominated by noise and so one usually wants to specify skip_levels to be 1 or 2. If skip_levels is 0 then all levels will be used to compute keypoint energy.

The upsample_highpasses and upsample_keypoint_energy parameters are used to control whether the individual subband coefficients and/org the keypoint energy map are upscaled by 2 before finding keypoints. If these parameters are None then no corresponding upscaling is performed. If non-None they specify the upscale method as outlined in dtcwt.sampling.upsample().

If method is None, the default 'fauqueur' method is used.

Name

Description

Parameters used

fauqueur

Geometric mean of absolute values[1]

alpha, beta

bendale

Minimum absolute value[2]

none

kingsbury

Cross-product of orthogonal highpasses

kappa

[1] Julien Fauqueur, Nick Kingsbury, and Ryan Anderson. Multiscale Keypoint Detection using the Dual-Tree Complex Wavelet Transform. 2006 International Conference on Image Processing, pages 1625-1628, October 2006. ISSN 1522-4880. doi: 10.1109/ICIP.2006.312656. http://ieeexplore.ieee.org/lpdocs/epic03/wrapper.htm?arnumber=4106857.

[2] Pashmina Bendale, Bill Triggs, and Nick Kingsbury. Multiscale Keypoint Analysis based on Complex Wavelets. In British Machine Vision Con-ference (BMVC), 2010. http://www-sigproc.eng.cam.ac.uk/~pb397/publications/BTK_BMVC_2010_abstract.pdf.

Image sampling

This module contains function for rescaling and re-sampling high- and low-pass highpasses.

Note

All of these functions take an integer co-ordinate (x, y) to be the centre of the corresponding pixel. Therefore the upper-left pixel notionally covers the interval (-0.5, 0.5) in x and y. An image with N rows and M columns, therefore, has an extent (-0.5, M-0.5) on the x-axis and an extent of (-0.5, N-0.5) on the y-axis. The rescale and upsample functions in this module will use this region as the extent of the image.

dtcwt.sampling.rescale(im, shape, method=None)

Return a resampled version of im scaled to shape.

Since the centre of pixel (x,y) has co-ordinate (x,y) the extent of im is actually \(x \in (-0.5, w-0.5]\) and \(y \in (-0.5, h-0.5]\) where (y,x) is im.shape. This returns a sampled version of im that has the same extent as a shape-sized array.

dtcwt.sampling.rescale_highpass(im, shape, method=None, sbs=None)

As rescale() except that the highpass image is first phase shifted to be centred on approximately DC, and has additional ‘sbs’ input allowing selection and re-ordering of subbands. This is useful mainly when the exact locations one wishes to sample from vary by subband.

‘sbs’ should ordinarily be a list of subband indices, in ascending order, e.g., np.array([0,2,3,5]) for just subbands 0, 2, 3 and 5; The returned array will be flattened to just 4 subbands. Pass [0,1,2,3,4,5] for all subbands.

Take care when re-ordering, preferably keeping the ‘sbs’ array outside the scope of this function to keep track of the new indices.

    1. Forshaw, Feb 2014.

dtcwt.sampling.sample(im, xs, ys, method=None)

Sample image at (x,y) given by elements of xs and ys. Both xs and ys must have identical shape and output will have this same shape. The location (x,y) refers to the centre of im[y,x]. Samples at fractional locations are calculated using the method specified by method (or 'lanczos' if method is None.)

Parameters
  • im – array to sample from

  • xs – x co-ordinates to sample

  • ys – y co-ordinates to sample

  • method – one of ‘bilinear’, ‘lanczos’ or ‘nearest’

Raises

ValueError – if xs and ys have differing shapes

dtcwt.sampling.sample_highpass(im, xs, ys, method=None, sbs=None)

As sample() except that the highpass image is first phase shifted to be centred on approximately DC, and has additional ‘sbs’ input allowing selection and re-ordering of subbands. This is useful mainly when the exact locations one wishes to sample from vary by subband.

‘sbs’ should ordinarily be a numpy array of subband indices, in ascending order, e.g., np.array([0,2,3,5]) for just subbands 0, 2, 3 and 5; The returned array will be flattened to just 4 subbands. Pass [0,1,2,3,4,5] for all subbands.

Take care when re-ordering, preferably keeping the ‘sbs’ array outside the scope of this function to keep track of the new indices.

    1. Forshaw, Feb 2014.

dtcwt.sampling.upsample(image, method=None)

Specialised function to upsample an image by a factor of two using a specified sampling method. If image is an array of shape (NxMx…) then the output will have shape (2Nx2Mx…). Only rows and columns are upsampled, depth axes and greater are interpolated but are not upsampled.

Parameters
  • image – an array containing the image to upsample

  • method – if non-None, a string specifying the sampling method to use.

If method is None, the default sampling method 'lanczos' is used. The following sampling methods are supported:

Name

Description

nearest

Nearest-neighbour sampling

bilinear

Bilinear sampling

lanczos

Lanczos sampling with window radius of 3

dtcwt.sampling.upsample_highpass(im, method=None)

As upsample() except that the highpass image is first phase rolled so that the filter has approximate DC centre frequency. The upshot is that this is the function to use when re-sampling complex subband images.

Image registration

Note

This module is experimental. It’s API may change between versions.

This module implements function for DTCWT-based image registration as outlined in [1]. These functions are 2D-only for the moment.

dtcwt.registration.estimatereg(source, reference, regshape=None, levels=None)

Estimate registration from which will map source to reference.

Parameters
  • source – transformed source image

  • reference – transformed reference image

The reference and source parameters should support the same API as dtcwt.Pyramid.

The local affine distortion is estimated at at 8x8 pixel scales. Return a NxMx6 array where the 6-element vector at (N,M) corresponds to the affine distortion parameters for the 8x8 block with index (N,M).

Use the velocityfield() function to convert the return value from this function into a velocity field.

If not-None, levels is a sequence of sequences of 0-based level indices to use when calculating the registration. If None then a default set of levels are used.

dtcwt.registration.velocityfield(avecs, shape, method=None)

Given the affine distortion parameters returned from estimatereg(), return a tuple of 2D arrays giving the x- and y- components of the velocity field. The shape of the velocity component field is shape. The velocities are measured in terms of normalised units where the image has width and height of unity.

The method parameter is interpreted as in dtcwt.sampling.rescale() and is the sampling method used to resize avecs to shape.

dtcwt.registration.warp(I, avecs, method=None)

A convenience function to warp an image according to the velocity field implied by avecs.

dtcwt.registration.warptransform(t, avecs, levels, method=None)

Return a warped version of a transformed image acting only on specified levels.

Parameters
  • t – a transformed image

  • avecs – an array of affine distortion parameters

  • levels – a sequence of 0-based indices specifying which levels to act on

t should be a dtcwt.Pyramid-compatible instance.

The method parameter is interpreted as in dtcwt.sampling.rescale() and is the sampling method used to resize avecs to shape.

Note

This function will clone the transform t but it is a shallow clone where possible. Only the levels specified in levels will be deep-copied and warped.

Plotting functions

Convenience functions for plotting DTCWT-related objects.

dtcwt.plotting.overlay_quiver(image, vectorField, level, offset)

Overlays nicely coloured quiver plot of complex coefficients over original full-size image, providing a useful phase visualisation.

Parameters
  • image – array holding grayscale values on the interval [0, 255] to display

  • vectorField – a single [MxNx6] numpy array of DTCWT coefficients

  • level – the transform level (1-indexed) of vectorField.

  • offset – Offset for DTCWT coefficients (typically 0.5)

Note

The level parameter is 1-indexed meaning that the third level has index “3”. This is unusual in Python but is kept for compatibility with similar MATLAB routines.

Should also work with other types of complex arrays (e.g., SLP coefficients), as long as the format is the same.

Usage example:

import dtcwt import dtcwt.plotting as plotting

mandrill = datasets.mandrill()

transform2d = dtcwt.Transform2d() mandrill_t = transform2d.forward(mandrill, nlevels=5)

plotting.overlay_quiver(mandrill*255, mandrill_t.highpasses[-1], 5, 0.5)

Miscellaneous and low-level support functions

Useful utilities for testing the 2-D DTCWT with synthetic images

dtcwt.utils.appropriate_complex_type_for(X)

Return an appropriate complex data type depending on the type of X. If X is already complex, return that, if it is floating point return a complex type of the appropriate size and if it is integer, choose an complex floating point type depending on the result of numpy.asfarray().

dtcwt.utils.as_column_vector(v)

Return v as a column vector with shape (N,1).

dtcwt.utils.asfarray(X)

Similar to numpy.asfarray() except that this function tries to preserve the original datatype of X if it is already a floating point type and will pass floating point arrays through directly without copying.

dtcwt.utils.drawcirc(r, w, du, dv, N)

Generate an image of size N*N pels, containing a circle radius r pels and centred at du,dv relative to the centre of the image. The edge of the circle is a cosine shaped edge of width w (from 10 to 90% points).

Python implementation by S. C. Forshaw, November 2013.

dtcwt.utils.drawedge(theta, r, w, N)

Generate an image of size N * N pels, of an edge going from 0 to 1 in height at theta degrees to the horizontal (top of image = 1 if angle = 0). r is a two-element vector, it is a coordinate in ij coords through which the step should pass. The shape of the intensity step is half a raised cosine w pels wide (w>=1).

T. E . Gale’s enhancement to drawedge() for MATLAB, transliterated to Python by S. C. Forshaw, Nov. 2013.

dtcwt.utils.reflect(x, minx, maxx)

Reflect the values in matrix x about the scalar values minx and maxx. Hence a vector x containing a long linearly increasing series is converted into a waveform which ramps linearly up and down between minx and maxx. If x contains integers and minx and maxx are (integers + 0.5), the ramps will have repeated max and min samples.

dtcwt.utils.stacked_2d_matrix_matrix_prod(mats1, mats2)

Interpret mats1 and mats2 as arrays of 2D matrices. I.e. mats1 has shape PxQxNxM and mats2 has shape PxQxMxR. The result is a PxQxNxR array equivalent to:

result[i,j,:,:] = mats1[i,j,:,:].dot(mats2[i,j,:,:])

for all valid row and column indices i and j.

dtcwt.utils.stacked_2d_matrix_vector_prod(mats, vecs)

Interpret mats and vecs as arrays of 2D matrices and vectors. I.e. mats has shape PxQxNxM and vecs has shape PxQxM. The result is a PxQxN array equivalent to:

result[i,j,:] = mats[i,j,:,:].dot(vecs[i,j,:])

for all valid row and column indices i and j.

dtcwt.utils.stacked_2d_vector_matrix_prod(vecs, mats)

Interpret mats and vecs as arrays of 2D matrices and vectors. I.e. mats has shape PxQxNxM and vecs has shape PxQxN. The result is a PxQxM array equivalent to:

result[i,j,:] = mats[i,j,:,:].T.dot(vecs[i,j,:])

for all valid row and column indices i and j.

dtcwt.utils.unpack(pyramid, backend='numpy')

Unpacks a pyramid give back the constituent parts.

Parameters
  • pyramid – The Pyramid of DTCWT transforms you wish to unpack

  • backend (str) – A string from ‘numpy’, ‘opencl’, or ‘tf’ indicating which attributes you want to unpack from the pyramid.

Returns

returns a generator which can be unpacked into the Yl, Yh and Yscale components of the pyramid. The generator will only return 2 values if the pyramid was created with the include_scale parameter set to false.

Note

You can still unpack a tf or opencl pyramid as if it were created by a numpy. In this case it will return a numpy array, rather than the backend specific array type.

Compatibility with MATLAB

Functions for compatibility with MATLAB scripts. These functions are intentionally similar in name and behaviour to the original functions from the DTCWT MATLAB toolbox. They are included in the library to ease the porting of MATLAB scripts but shouldn’t be used in new projects.

Note

The functionality of dtwavexfm2b and dtwaveifm2b has been folded into dtwavexfm2 and dtwaveifm2. For convenience of porting MATLAB scripts, the original function names are available in the dtcwt module as aliases but they should not be used in new code.

dtcwt.compat.dtwaveifm(Yl, Yh, biort='near_sym_a', qshift='qshift_a', gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 1D reconstruction.

Parameters
  • Yl – The real lowpass subband from the final level

  • Yh – A sequence containing the complex highpass subband for each level.

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • gain_mask – Gain to be applied to each subband.

Returns Z

Reconstructed real array.

The l-th element of gain_mask is gain for wavelet subband at level l. If gain_mask[l] == 0, no computation is performed for band l. Default gain_mask is all ones. Note that l is 0-indexed.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Example:

# Performs a reconstruction from Yl,Yh using the 13,19-tap filters
# for level 1 and the Q-shift 14-tap filters for levels >= 2.
Z = dtwaveifm(Yl, Yh, 'near_sym_b', 'qshift_b')
dtcwt.compat.dtwaveifm2(Yl, Yh, biort='near_sym_a', qshift='qshift_a', gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 2D reconstruction.

Parameters
  • Yl – The real lowpass subband from the final level

  • Yh – A sequence containing the complex highpass subband for each level.

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • gain_mask – Gain to be applied to each subband.

Returns Z

Reconstructed real array

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Example:

# Performs a 3-level reconstruction from Yl,Yh using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Z = dtwaveifm2(Yl, Yh, 'near_sym_b', 'qshift_b')
dtcwt.compat.dtwaveifm2b(Yl, Yh, biort='near_sym_a', qshift='qshift_a', gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 2D reconstruction.

Parameters
  • Yl – The real lowpass subband from the final level

  • Yh – A sequence containing the complex highpass subband for each level.

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • gain_mask – Gain to be applied to each subband.

Returns Z

Reconstructed real array

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Example:

# Performs a 3-level reconstruction from Yl,Yh using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Z = dtwaveifm2(Yl, Yh, 'near_sym_b', 'qshift_b')
dtcwt.compat.dtwaveifm3(Yl, Yh, biort='near_sym_a', qshift='qshift_a', ext_mode=4)

Perform an n-level dual-tree complex wavelet (DTCWT) 3D reconstruction.

Parameters
  • Yl – The real lowpass subband from the final level

  • Yh – A sequence containing the complex highpass subband for each level.

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • ext_mode – Extension mode. See below.

Returns Z

Reconstructed real image matrix.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

There are two values for ext_mode, either 4 or 8. If ext_mode = 4, check whether 1st level is divisible by 2 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coefs can be divided by 4. If any dimension size is not a multiple of 4, append extra coefs by repeating the edges. If ext_mode = 8, check whether 1st level is divisible by 4 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coeffs can be divided by 8. If any dimension size is not a multiple of 8, append extra coeffs by repeating the edges twice.

Example:

# Performs a 3-level reconstruction from Yl,Yh using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Z = dtwaveifm3(Yl, Yh, 'near_sym_b', 'qshift_b')
dtcwt.compat.dtwavexfm(X, nlevels=3, biort='near_sym_a', qshift='qshift_a', include_scale=False)

Perform a n-level DTCWT decompostion on a 1D column vector X (or on the columns of a matrix X).

Parameters
  • X – 1D real array or 2D real array whose columns are to be transformed

  • nlevels – Number of levels of wavelet decomposition

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

Returns Yl

The real lowpass image from the final level

Returns Yh

A tuple containing the (N, M, 6) shape complex highpass subimages for each level.

Returns Yscale

If include_scale is True, a tuple containing real lowpass coefficients for every scale.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Example:

# Performs a 5-level transform on the real image X using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Yl, Yh = dtwavexfm(X,5,'near_sym_b','qshift_b')
dtcwt.compat.dtwavexfm2(X, nlevels=3, biort='near_sym_a', qshift='qshift_a', include_scale=False)

Perform a n-level DTCWT-2D decompostion on a 2D matrix X.

Parameters
Returns Yl

The real lowpass image from the final level

Returns Yh

A tuple containing the complex highpass subimages for each level.

Returns Yscale

If include_scale is True, a tuple containing real lowpass coefficients for every scale.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Example:

# Performs a 3-level transform on the real image X using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Yl, Yh = dtwavexfm2(X, 3, 'near_sym_b', 'qshift_b')
dtcwt.compat.dtwavexfm2b(X, nlevels=3, biort='near_sym_a', qshift='qshift_a', include_scale=False)

Perform a n-level DTCWT-2D decompostion on a 2D matrix X.

Parameters
Returns Yl

The real lowpass image from the final level

Returns Yh

A tuple containing the complex highpass subimages for each level.

Returns Yscale

If include_scale is True, a tuple containing real lowpass coefficients for every scale.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Example:

# Performs a 3-level transform on the real image X using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Yl, Yh = dtwavexfm2(X, 3, 'near_sym_b', 'qshift_b')
dtcwt.compat.dtwavexfm3(X, nlevels=3, biort='near_sym_a', qshift='qshift_a', include_scale=False, ext_mode=4, discard_level_1=False)

Perform a n-level DTCWT-3D decompostion on a 3D matrix X.

Parameters
  • X – 3D real array-like object

  • nlevels – Number of levels of wavelet decomposition

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • ext_mode – Extension mode. See below.

  • discard_level_1 – True if level 1 high-pass bands are to be discarded.

Returns Yl

The real lowpass image from the final level

Returns Yh

A tuple containing the complex highpass subimages for each level.

Each element of Yh is a 4D complex array with the 4th dimension having size 28. The 3D slice Yh[l][:,:,:,d] corresponds to the complex higpass coefficients for direction d at level l where d and l are both 0-indexed.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

There are two values for ext_mode, either 4 or 8. If ext_mode = 4, check whether 1st level is divisible by 2 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coefs can be divided by 4. If any dimension size is not a multiple of 4, append extra coefs by repeating the edges. If ext_mode = 8, check whether 1st level is divisible by 4 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coeffs can be divided by 8. If any dimension size is not a multiple of 8, append extra coeffs by repeating the edges twice.

If discard_level_1 is True the highpass coefficients at level 1 will be discarded. (And, in fact, will never be calculated.) This turns the transform from being 8:1 redundant to being 1:1 redundant at the cost of no-longer allowing perfect reconstruction. If this option is selected then Yh[0] will be None. Note that dtwaveifm3() will accepts Yh[0] being None and will treat it as being zero.

Example:

# Performs a 3-level transform on the real 3D array X using the 13,19-tap
# filters for level 1 and the Q-shift 14-tap filters for levels >= 2.
Yl, Yh = dtwavexfm3(X, 3, 'near_sym_b', 'qshift_b')

Backends

The following modules provide backend-specific implementations. Usually you won’t need to import these modules directly; the main API will use an appropriate implementation. Occasionally, however, you may want to benchmark one implementation against the other.

NumPy

A backend which uses NumPy to perform the filtering. This backend should always be available.

class dtcwt.numpy.Pyramid(lowpass, highpasses, scales=None)

A representation of a transform domain signal.

Backends are free to implement any class which respects this interface for storing transform-domain signals. The inverse transform may accept a backend-specific version of this class but should always accept any class which corresponds to this interface.

lowpass

A NumPy-compatible array containing the coarsest scale lowpass signal.

highpasses

A tuple where each element is the complex subband coefficients for corresponding scales finest to coarsest.

scales

(optional) A tuple where each element is a NumPy-compatible array containing the lowpass signal for corresponding scales finest to coarsest. This is not required for the inverse and may be None.

class dtcwt.numpy.Transform1d(biort='near_sym_a', qshift='qshift_a')

An implementation of the 1D DT-CWT in NumPy.

Parameters
forward(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT decompostion on a 1D column vector X (or on the columns of a matrix X).

Parameters
  • X – 1D real array or 2D real array whose columns are to be transformed

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.Pyramid-like object representing the transform result.

If biort or qshift are strings, they are used as an argument to the biort() or qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

inverse(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 1D reconstruction.

Parameters
  • pyramid – A dtcwt.Pyramid-like object containing the transformed signal.

  • gain_mask – Gain to be applied to each subband.

Returns

Reconstructed real array.

The l-th element of gain_mask is gain for wavelet subband at level l. If gain_mask[l] == 0, no computation is performed for band l. Default gain_mask is all ones. Note that l is 0-indexed.

class dtcwt.numpy.Transform2d(biort='near_sym_a', qshift='qshift_a')

An implementation of the 2D DT-CWT via NumPy. biort and qshift are the wavelets which parameterise the transform.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

forward(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT-2D decompostion on a 2D matrix X.

Parameters
  • X – 2D real array

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.Pyramid compatible object representing the transform-domain signal

inverse(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 2D reconstruction.

Parameters
  • pyramid – A dtcwt.Pyramid-like class holding the transform domain representation to invert.

  • gain_mask – Gain to be applied to each subband.

Returns

A numpy-array compatible instance with the reconstruction.

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

class dtcwt.numpy.Transform3d(biort='near_sym_a', qshift='qshift_a', ext_mode=4)

An implementation of the 3D DT-CWT via NumPy. biort and qshift are the wavelets which parameterise the transform. Valid values are documented in dtcwt.coeffs.biort() and dtcwt.coeffs.qshift().

forward(X, nlevels=3, include_scale=False, discard_level_1=False)

Perform a n-level DTCWT-3D decompostion on a 3D matrix X.

Parameters
  • X – 3D real array-like object

  • nlevels – Number of levels of wavelet decomposition

  • biort – Level 1 wavelets to use. See dtcwt.coeffs.biort().

  • qshift – Level >= 2 wavelets to use. See dtcwt.coeffs.qshift().

  • discard_level_1 – True if level 1 high-pass bands are to be discarded.

Returns

a dtcwt.Pyramid instance

Each element of the Pyramid highpasses tuple is a 4D complex array with the 4th dimension having size 28. The 3D slice [l][:,:,:,d] corresponds to the complex higpass coefficients for direction d at level l where d and l are both 0-indexed.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

There are two values for ext_mode, either 4 or 8. If ext_mode = 4, check whether 1st level is divisible by 2 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coefs can be divided by 4. If any dimension size is not a multiple of 4, append extra coefs by repeating the edges. If ext_mode = 8, check whether 1st level is divisible by 4 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coeffs can be divided by 8. If any dimension size is not a multiple of 8, append extra coeffs by repeating the edges twice.

If discard_level_1 is True the highpass coefficients at level 1 will not be discarded. (And, in fact, will never be calculated.) This turns the transform from being 8:1 redundant to being 1:1 redundant at the cost of no-longer allowing perfect reconstruction. If this option is selected then the first element of the highpasses tuple will be None. Note that dtcwt.Transform3d.inverse() will accept the first element being None and will treat it as being zero.

inverse(pyramid)

Perform an n-level dual-tree complex wavelet (DTCWT) 3D reconstruction.

Parameters
  • pyramid – The dtcwt.Pyramid-like instance representing the transformed signal.

  • biort – Level 1 wavelets to use. See biort().

  • qshift – Level >= 2 wavelets to use. See qshift().

  • ext_mode – Extension mode. See below.

Returns

Reconstructed real image matrix.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

There are two values for ext_mode, either 4 or 8. If ext_mode = 4, check whether 1st level is divisible by 2 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coefs can be divided by 4. If any dimension size is not a multiple of 4, append extra coefs by repeating the edges. If ext_mode = 8, check whether 1st level is divisible by 4 (if not we raise a ValueError). Also check whether from 2nd level onwards, the coeffs can be divided by 8. If any dimension size is not a multiple of 8, append extra coeffs by repeating the edges twice.

dtcwt.numpy.lowlevel.coldfilt(X, ha, hb)

Filter the columns of image X using the two filters ha and hb = reverse(ha). ha operates on the odd samples of X and hb on the even samples. Both filters should be even length, and h should be approx linear phase with a quarter sample advance from its mid pt (i.e. \(|h(m/2)| > |h(m/2 + 1)|\)).

                  ext        top edge                     bottom edge       ext
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a

The output is decimated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X columns before each filter is applied.

Raises ValueError if the number of rows in X is not a multiple of 4, the length of ha does not match hb or the lengths of ha or hb are non-even.

dtcwt.numpy.lowlevel.colfilter(X, h)

Filter the columns of image X using filter vector h, without decimation. If len(h) is odd, each output sample is aligned with each input sample and Y is the same size as X. If len(h) is even, each output sample is aligned with the mid point of each pair of input samples, and Y.shape = X.shape + [1 0].

Parameters
  • X – an image whose columns are to be filtered

  • h – the filter coefficients.

Returns Y

the filtered image.

dtcwt.numpy.lowlevel.colifilt(X, ha, hb)

Filter the columns of image X using the two filters ha and hb = reverse(ha). ha operates on the odd samples of X and hb on the even samples. Both filters should be even length, and h should be approx linear phase with a quarter sample advance from its mid pt (i.e :math:`|h(m/2)| > |h(m/2 + 1)|).

                  ext       left edge                      right edge       ext
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a

The output is interpolated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X columns before each filter is applied.

OpenCL

Provide low-level OpenCL accelerated operations. This backend requires that PyOpenCL be installed.

class dtcwt.opencl.Pyramid(lowpass, highpasses, scales=None)

An interface-compatible version of dtcwt.Pyramid where the initialiser arguments are assumed to by pyopencl.array.Array instances.

The attributes defined in dtcwt.Pyramid are implemented via properties. The original OpenCL arrays may be accessed via the cl_... attributes.

Note

The copy from device to host is performed once and then memoized. This makes repeated access to the host-side attributes efficient but will mean that any changes to the device-side arrays will not be reflected in the host-side attributes after their first access. You should not be modifying the arrays once you return an instance of this class anyway but if you do, beware!

cl_lowpass

The CL array containing the lowpass image.

cl_highpasses

A tuple of CL arrays containing the subband images.

cl_scales

(optional) Either None or a tuple of lowpass images for each scale.

class dtcwt.opencl.Transform2d(biort='near_sym_a', qshift='qshift_a', queue=None)

An implementation of the 2D DT-CWT via OpenCL. biort and qshift are the wavelets which parameterise the transform.

If queue is non-None it is an instance of pyopencl.CommandQueue which is used to compile and execute the OpenCL kernels which implement the transform. If it is None, the first available compute device is used.

If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Note

At the moment only the forward transform is accelerated. The inverse transform uses the NumPy backend.

forward(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT-2D decompostion on a 2D matrix X.

Parameters
  • X – 2D real array

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.Pyramid compatible object representing the transform-domain signal

Note

X may be a pyopencl.array.Array instance which has already been copied to the device. In which case, it must be 2D. (I.e. a vector will not be auto-promoted.)

inverse(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 2D reconstruction.

Parameters
  • pyramid – A dtcwt.Pyramid-like class holding the transform domain representation to invert.

  • gain_mask – Gain to be applied to each subband.

Returns

A numpy-array compatible instance with the reconstruction.

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

exception dtcwt.opencl.lowlevel.NoCLPresentError
dtcwt.opencl.lowlevel.axis_convolve(X, h, axis=0, queue=None, output=None)

Filter along an of X using filter vector h. If h has odd length, each output sample is aligned with each input sample and Y is the same size as X. If h has even length, each output sample is aligned with the mid point of each pair of input samples, and the output matrix’s shape is increased by one along the convolution axis.

After convolution, the pyopencl.array.Array instance holding the device-side output is returned. This may be accessed on the host via to_array().

The axis of convolution is specified by axis. The default direction of convolution is column-wise.

If queue is non-None, it should be a pyopencl.CommandQueue instance which is used to perform the computation. If None, a default global queue is used.

If output is non-None, it should be a pyopencl.array.Array instance which the result is written into. If None, an output array is created.

dtcwt.opencl.lowlevel.coldfilt(X, ha, hb, queue=None)

Filter the columns of image X using the two filters ha and hb = reverse(ha). ha operates on the odd samples of X and hb on the even samples. Both filters should be even length, and h should be approx linear phase with a quarter sample advance from its mid pt (i.e. \(|h(m/2)| > |h(m/2 + 1)|\)).

                  ext        top edge                     bottom edge       ext
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a

The output is decimated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X columns before each filter is applied.

Raises ValueError if the number of rows in X is not a multiple of 4, the length of ha does not match hb or the lengths of ha or hb are non-even.

dtcwt.opencl.lowlevel.colfilter(X, h)

Filter the columns of image X using filter vector h, without decimation. If len(h) is odd, each output sample is aligned with each input sample and Y is the same size as X. If len(h) is even, each output sample is aligned with the mid point of each pair of input samples, and Y.shape = X.shape + [1 0].

The filtering will be accelerated via OpenCL.

Parameters
  • X – an image whose columns are to be filtered

  • h – the filter coefficients.

Returns Y

the filtered image.

dtcwt.opencl.lowlevel.colifilt(X, ha, hb, queue=None)

Filter the columns of image X using the two filters ha and hb = reverse(ha). ha operates on the odd samples of X and hb on the even samples. Both filters should be even length, and h should be approx linear phase with a quarter sample advance from its mid pt (i.e :math:`|h(m/2)| > |h(m/2 + 1)|).

                  ext       left edge                      right edge       ext
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a

The output is interpolated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X columns before each filter is applied.

dtcwt.opencl.lowlevel.get_default_queue()

Return the default queue used for computation if one is not specified.

This function is memoized and so only one queue is created after multiple invocations.

Tensorflow

Currently the Tensorflow backend only supports single precision operations, and only has functionality for the Transform1d() and Transform2d() classes (i.e. changing the backend to ‘tf’ will still use the numpy Transform3d() class).

To preserve functionality, the Transform1d() and Transform2d() classes have a forward method which behaves identically to the NumPy backend. However, to get speedups with tensorflow, we want to feed our transform batches of images. For this reason, the 1-D and 2-D transforms also have forward_channels and inverse_channels methods. See the below documentation for how to use these.

Provide low-level Tensorflow accelerated operations. This backend requires that Tensorflow be installed. Works best with a GPU but still offers good improvements with a CPU.

class dtcwt.tf.Pyramid(lowpass, highpasses, scales=None, numpy=False)

A tensorflow representation of a transform domain signal.

An interface-compatible version of dtcwt.Pyramid where the initialiser arguments are assumed to be tf.Variable instances.

The attributes defined in dtcwt.Pyramid are implemented via properties. The original tf arrays may be accessed via the ..._op(s) attributes.

lowpass_op

A tensorflow tensor that can be evaluated in a session to return the coarsest scale lowpass signal for the input, X.

highpasses_op

A tuple of tensorflow tensors, where each element is the complex subband coefficients for corresponding scales finest to coarsest.

scales_ops

(optional) A tuple where each element is a tensorflow tensor containing the lowpass signal for corresponding scales finest to coarsest. This is not required for the inverse and may be None.

class dtcwt.tf.Transform1d(biort='near_sym_a', qshift='qshift_a')

An implementation of the 1D DT-CWT in Tensorflow.

Parameters

Note

Calling the methods in this class with different inputs will slightly vary the results. If you call the forward() or forward_channels() methods with a numpy array, they load this array into a tf.Variable and create the graph. Subsequent calls to dtcwt.tf.Pyramid.lowpass or other attributes in the pyramid will create a session and evaluate these parameters. If the above methods are called with a tensorflow variable or placeholder, these will be used to create the graph. As such, to evaluate the results, you will need to look at the dtcwt.tf.Pyramid.lowpass_op attribute (calling the lowpass attribute will try to evaluate the graph with no initialized variables and likely result in a runtime error).

The behaviour is similar for the inverse() and inverse_channels() methods, except these return an array, rather than a Pyramid style class. If a dtcwt.tf.Pyramid was created by calling the forward methods with a numpy array, providing this pyramid to the inverse methods will return a numpy array. If however a dtcwt.tf.Pyramid was created by calling the forward methods with a tensorflow variable, the result from calling the inverse methods will also be a tensorflow variable.

forward(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT decompostion on a 1D column vector X (or on the columns of a matrix X).

Can provide the forward transform with either an np array (naive usage), or a tensorflow variable or placeholder (designed usage). To transform batches of vectors, use the forward_channels() method.

Parameters
  • X – 1D real array or 2D real array whose columns are to be transformed.

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.tf.Pyramid object representing the transform result.

If biort or qshift are strings, they are used as an argument to the biort() or qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

forward_channels(X, nlevels=3, include_scale=False)

Perform a n-level DTCWT decompostion on a 3D array X.

Can provide the forward transform with either an np array (naive usage), or a tensorflow variable or placeholder (designed usage).

Parameters
  • X – 3D real array. Batch of matrices whose columns are to be transformed (i.e. the second dimension).

  • nlevels – Number of levels of wavelet decomposition

Returns

A dtcwt.tf.Pyramid object representing the transform result.

If biort or qshift are strings, they are used as an argument to the biort() or qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

inverse(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 1D reconstruction.

Parameters
  • pyramid – A dtcwt.Pyramid-like object containing the transformed signal.

  • gain_mask – Gain to be applied to each subband.

Returns

Reconstructed real array. Will be a tf Variable if the Pyramid was made with tf inputs, otherwise a numpy array.

The l-th element of gain_mask is gain for wavelet subband at level l. If gain_mask[l] == 0, no computation is performed for band l. Default gain_mask is all ones. Note that l is 0-indexed.

inverse_channels(pyramid, gain_mask=None)

Perform an n-level dual-tree complex wavelet (DTCWT) 1D reconstruction on a 3D array of signals. The inverse is done on the second dimension of these.

This is designed to work after calling the forward_channels() method.

Parameters
  • pyramid – A dtcwt.Pyramid-like object containing the transformed signal. The lowpass signal in the pyramid should be a 3D array to use this method.

  • gain_mask – Gain to be applied to each subband.

Returns

Reconstructed array. Will be a tf Variable if the Pyramid was made with tf inputs, otherwise a numpy array.

The l-th element of gain_mask is gain for wavelet subband at level l. If gain_mask[l] == 0, no computation is performed for band l. Default gain_mask is all ones. Note that l is 0-indexed.

class dtcwt.tf.Transform2d(biort='near_sym_a', qshift='qshift_a')

An implementation of the 2D DT-CWT via Tensorflow.

Parameters
  • biort – The biorthogonal wavelet family to use.

  • qshift – The quarter shift wavelet family to use.

Note

biort and qshift are the wavelets which parameterise the transform. If biort or qshift are strings, they are used as an argument to the dtcwt.coeffs.biort() or dtcwt.coeffs.qshift() functions. Otherwise, they are interpreted as tuples of vectors giving filter coefficients. In the biort case, this should be (h0o, g0o, h1o, g1o). In the qshift case, this should be (h0a, h0b, g0a, g0b, h1a, h1b, g1a, g1b).

Note

Calling the methods in this class with different inputs will slightly vary the results. If you call the forward() or forward_channels() methods with a numpy array, they load this array into a tf.Variable and create the graph. Subsequent calls to dtcwt.tf.Pyramid.lowpass or other attributes in the pyramid will create a session and evaluate these parameters. If the above methods are called with a tensorflow variable or placeholder, these will be used to create the graph. As such, to evaluate the results, you will need to look at the dtcwt.tf.Pyramid.lowpass_op attribute (calling the lowpass attribute will try to evaluate the graph with no initialized variables and likely result in a runtime error).

The behaviour is similar for the inverse methods, except these return an array, rather than a Pyramid style class. If a dtcwt.tf.Pyramid was created by calling the forward methods with a numpy array, providing this pyramid to the inverse methods will return a numpy array. If however a dtcwt.tf.Pyramid was created by calling the forward methods with a tensorflow variable, the result from calling the inverse methods will also be a tensorflow variable.

forward(X, nlevels=3, include_scale=False)

Perform a forward transform on an image.

Can provide the forward transform with either an np array (naive usage), or a tensorflow variable or placeholder (designed usage). To transform batches of images, use the forward_channels() method.

Parameters
  • X (ndarray) – Input image which you wish to transform. Can be a numpy array, tensorflow Variable or tensorflow placeholder. See comments below.

  • nlevels (int) – Number of levels of the dtcwt transform to calculate.

  • include_scale (bool) – Whether or not to return the lowpass results at each scale of the transform, or only at the highest scale (as is custom for multi-resolution analysis)

Returns

A dtcwt.tf.Pyramid object

Note

If a numpy array is provided, the forward function will create a tensorflow variable to hold the input image, and then create the graph of the right size to match the input, and then feed the input into the graph and evaluate it. This operation will return a Pyramid object similar to how running the numpy version would.

forward_channels(X, data_format, nlevels=3, include_scale=False)

Perform a forward transform on an image with multiple channels.

Will perform the DTCWT independently on each channel.

Parameters
  • X – Input image which you wish to transform.

  • nlevels (int) – Number of levels of the dtcwt transform to calculate.

  • include_scale (bool) – Whether or not to return the lowpass results at each scale of the transform, or only at the highest scale (as is custom for multiresolution analysis)

  • data_format (str) –

    An optional string of the form: “nhw” (or “chw”), “hwn” (or “hwc”), “nchw” or “nhwc”. Note that for these strings, ‘n’ is used to indicate where the batch dimension is, ‘c’ is used to indicate where the image channels are, ‘h’ is used to indicate where the row dimension is, and ‘c’ is used to indicate where the columns are. If the data_format is:

    • ”nhw” : the input will be interpreted as a batch of 2D images, with the batch dimension as the first.

    • ”chw” : will function exactly the same as “nhw” but is offered to indicate the input is a 2D image with channels.

    • ”hwn” : the input will be interpreted as a batch of 2D images with the batch dimension as the last.

    • ”hwc” : will function exatly the same as “hwc” but is offered to indicate the input is a 2D image with channels.

    • ”nchw” : the input is a batch of images with channel dimension as the second dimension. Batch dimension is first.

    • ”nhwc” : the input is a batch of images with channel dimension as the last dimension. Batch dimension is first.

Returns

A dtcwt.tf.Pyramid object

inverse(pyramid, gain_mask=None)

Perform an inverse transform on an image.

Can provide the inverse transform with either an np array (naive usage), or a tensorflow variable or placeholder (designed usage).

Parameters
  • pyramid – A dtcwt.tf.Pyramid like class holding the transform domain representation to invert

  • gain_mask – Gain to be applied to each sub-band. Should have shape (6, nlevels) or be None.

Returns

An array , X, compatible with the reconstruction. Will be a tf Variable if the Pyramid was made with tf inputs, otherwise a numpy array.

Note

A tf.Variable is returned if the pyramid input was a Pyramid class. If it wasn’t, then, we return a numpy array (note that this is inefficient, as in both cases we have to construct the graph - in the second case, we then execute it and discard it).

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

inverse_channels(pyramid, data_format, gain_mask=None)

Perform an inverse transform on an image with multiple channels.

Must provide with a tensorflow variable or placeholder (unlike the more general inverse()).

This is designed to work after calling the forward_channels() method. You must use the same data_format for the inverse_channels as the one used for the forward_channels (unless you have explicitly reshaped the output).

Parameters
  • pyramid – A dtcwt.tf.Pyramid like class holding the transform domain representation to invert

  • data_format (str) –

    An optional string of the form: “nhw” (or “chw”), “hwn” (or “hwc”), “nchw” or “nhwc”. Note that for these strings, ‘n’ is used to indicate where the batch dimension is, ‘c’ is used to indicate where the image channels are, ‘h’ is used to indicate where the row dimension is, and ‘c’ is used to indicate where the columns are. If the data_format is:

    * "nhw" - the input will be interpreted as a batch of 2D images,
      with the batch dimension as the first.
    * "chw" - will function exactly the same as "nhw" but it offered
      to indicate the input is a 2D image with channels.
    * "hwn" - the input will be interpreted as a batch of 2D images
      with the batch dimension as the last.
    * "hwc" - will function exatly the same as "hwc" but is offered
      to indicate the input is a 2D image with channels.
    * "nchw" - the input is a batch of images with channel dimension
      as the second dimension. Batch dimension is first.
    * "nhwc" - the input is a batch of images with channel dimension
      as the last dimension. Batch dimension is first.
    

  • gain_mask – Gain to be applied to each subband. Should have shape [6, nlevels].

Returns

An array , X, compatible with the reconstruction. Will be a tf Variable if the Pyramid was made with tf inputs, otherwise a numpy array.

The (d, l)-th element of gain_mask is gain for subband with direction d at level l. If gain_mask[d,l] == 0, no computation is performed for band (d,l). Default gain_mask is all ones. Note that both d and l are zero-indexed.

dtcwt.tf.lowlevel.coldfilt(X, ha, hb, no_decimate=False)

Filter the columns of image X using the two filters ha and hb = reverse(ha).

Parameters
  • X – The input, of size [batch, h, w]

  • ha – Filter to be used on the odd samples of x.

  • hb – Filter to bue used on the even samples of x.

  • no_decimate – If true, keep the same input size

Both filters should be even length, and h should be approx linear phase with a quarter sample (i.e. an \(e^{j \pi/4}\)) advance from its mid pt (i.e. \(|h(m/2)| > |h(m/2 + 1)|\)):

                  ext        top edge                     bottom edge       ext
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a

The output is decimated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X columns before each filter is applied.

:raises ValueError if the number of rows in X is not a multiple of 4, the

length of ha does not match hb or the lengths of ha or hb are non-even.

dtcwt.tf.lowlevel.colfilter(X, h, align=False)

Filter the columns of image X using filter vector h, without decimation.

Parameters
  • X – an image whose columns are to be filtered

  • h – the filter coefficients.

  • align – If true, then will have Y keep the same output shape as X, even if h has even length. Makes no difference if len(h) is odd.

Returns Y

the filtered image.

If len(h) is odd, each output sample is aligned with each input sample and Y is the same size as X. If len(h) is even, each output sample is aligned with the mid point of each pair of input samples, and Y.shape = X.shape + [1 0].

dtcwt.tf.lowlevel.colifilt(X, ha, hb, no_decimate=False)

Filter the columns of image X using the two filters ha and hb = reverse(ha).

Parameters
  • X – The input, of size [batch, h, w]

  • ha – Filter to be used on the odd samples of x.

  • hb – Filter to bue used on the even samples of x.

  • no_decimate – Not implemented yet

Both filters should be even length, and h should be approx linear phase with a quarter sample advance from its mid pt (i.e :math:`|h(m/2)| > |h(m/2 + 1)|).

                  ext       left edge                      right edge       ext
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a

The output is interpolated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X columns before each filter is applied.

dtcwt.tf.lowlevel.rowdfilt(X, ha, hb, no_decimate=False)

Filter the rows of image X using the two filters ha and hb = reverse(ha).

Parameters
  • X – The input, of size [batch, h, w]

  • ha – Filter to be used on the odd samples of x.

  • hb – Filter to bue used on the even samples of x.

  • no_decimate – If true, keep the same input size

Both filters should be even length, and h should be approx linear phase with a quarter sample advance from its mid pt (i.e. \(|h(m/2)| > |h(m/2 + 1)|\)):

                  ext        top edge                     bottom edge       ext
Level 1:        !               |               !               |               !
odd filt on .    b   b   b   b   a   a   a   a   a   a   a   a   b   b   b   b
odd filt on .      a   a   a   a   b   b   b   b   b   b   b   b   a   a   a   a
Level 2:        !               |               !               |               !
+q filt on x      b       b       a       a       a       a       b       b
-q filt on o          a       a       b       b       b       b       a       a

The output is decimated by two from the input sample rate and the results from the two filters, Ya and Yb, are interleaved to give Y. Symmetric extension with repeated end samples is used on the composite X rows before each filter is applied.

:raises ValueError if the number of columns in X is not a multiple of 4, the

length of ha does not match hb or the lengths of ha or hb are non-even.

dtcwt.tf.lowlevel.rowfilter(X, h, align=False)

Filter the rows of image X using filter vector h, without decimation.

Parameters
  • X – a tensor of images whose rows are to be filtered

  • h – the filter coefficients.

  • align – If true, then will have Y keep the same output shape as X, even if h has even length. Makes no difference if len(h) is odd.

Returns Y

the filtered image.

If len(h) is odd, each output sample is aligned with each input sample and Y is the same size as X. If len(h) is even, each output sample is aligned with the mid point of each pair of input samples, and Y.shape = X.shape + [0 1].