API Reference

Computing the DT-CWT

dtcwt.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 biort().
  • qshift – Level >= 2 wavelets to use. See 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 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).

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.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 biort().
  • qshift – Level >= 2 wavelets to use. See 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 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).

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.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:
  • X – 2D real array
  • nlevels – Number of levels of wavelet decomposition
  • biort – Level 1 wavelets to use. See biort().
  • qshift – Level >= 2 wavelets to use. See qshift().
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 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).

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.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 biort().
  • qshift – Level >= 2 wavelets to use. See 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 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).

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.dtwavexfm3(X, nlevels=3, biort='near_sym_a', qshift='qshift_a', 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 biort().
  • qshift – Level >= 2 wavelets to use. See 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 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).

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')
dtcwt.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 biort().
  • qshift – Level >= 2 wavelets to use. See 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 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).

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

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

Raises:
  • IOError – if name does not correspond to a set of wavelets known to the library.
  • ValueError – if name specifies a qshift() wavelet.
dtcwt.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.

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

Raises:
  • IOError – if name does not correspond to a set of wavelets known to the library.
  • ValueError – if name specifies a biort() wavelet.

Keypoint analysis

dtcwt.keypoint.find_keypoints(highpass_subbands, method=None, alpha=1.0, beta=0.4, kappa=0.16666666666666666, threshold=None, max_points=None, upsample_keypoint_energy=None, upsample_subbands=None, refine_positions=True, skip_levels=1)
Parameters:
  • highpass_subbands – (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_subbands 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 subbands 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

Rescaling and re-sampling high- and low-pass subbands.

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)

As sample() except that the highpass image is first phase shifted to be centred on approximately DC.

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)

As rescale() except that the highpass image is first phase shifted to be centred on approximately DC.

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.

Low-level support functions

A normal user should not need to call these functions but they are documented here just in case you do.

dtcwt.lowlevel.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.lowlevel.as_column_vector(v)

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

dtcwt.lowlevel.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.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.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.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.

dtcwt.lowlevel.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.

Table Of Contents

Previous topic

Example scripts

This Page