API#

curvelops#

curvelops.curvelops#

Provides a LinearOperator for the 2D and 3D curvelet transforms.

class curvelops.curvelops.FDCT(dims, axes, nbscales=None, nbangles_coarse=16, allcurvelets=True, dtype='complex128')[source]#

Bases: LinearOperator

2D/3D dimensional Curvelet operator. Apply 2D/3D Curvelet Transform along two axes of a multi-dimensional array of size dims.

Parameters#

dimstuple

Number of samples for each dimension.

axestuple, optional

Axes along which FDCT is applied.

nbscalesint, optional

Number of scales (including the coarsest level); Defaults to ceil(log2(min(input_dims)) - 3).

nbangles_coarseint, optional

Number of angles at 2nd coarsest scale.

allcurveletsbool, optional

Use curvelets at the finest (last) scale. If False, a wavelet transform will be used for the finest scale. The coarsest scale is always a wavelet transform; the ones between the coarsest and the finest are all curvelet transforms. This option only affects the finest scale.

dtypeDTypeLike, optional

dtype of the transform.

inverse(x)[source]#
struct(x)[source]#

Convert curvelet flattened vector to curvelet structure.

The FDCT always returns a 1D vector that has all curvelet coefficients. These coefficients can be organized into scales, wedges and spatial positions. Applying this function to a 1D vector generates this structure.

Parameters#

xNDArray

Input flattened vector.

Returns#

FDCTStructLike

Curvelet structure, a list of lists of multidimensional arrays. The first index corresponds to scale, the second corresponds to angular wedge.

vect(x)[source]#

Convert curvelet structure to curvelet flattened vector.

The FDCT always returns a 1D vector that has all curvelet coefficients. These coefficients can be organized into scales, wedges and spatial positions. Applying this function to a curvelet structure returns the flattened vector.

Parameters#

xFDCTStructLike

Input curvelet structure.

Returns#

NDArray

Flattened vector.

class curvelops.curvelops.FDCT2D(dims, axes=(-2, -1), nbscales=None, nbangles_coarse=16, allcurvelets=True, dtype='complex128')[source]#

Bases: FDCT

2D dimensional Curvelet operator. Apply 2D Curvelet Transform along two axes of a multi-dimensional array of size dims.

Parameters#

dimstuple

Number of samples for each dimension.

axestuple, optional

Axes along which FDCT is applied.

nbscalesint, optional

Number of scales (including the coarsest level); Defaults to ceil(log2(min(input_dims)) - 3).

nbangles_coarseint, optional

Number of angles at 2nd coarsest scale.

allcurveletsbool, optional

Use curvelets at the finest (last) scale. If False, a wavelet transform will be used for the finest scale. The coarsest scale is always a wavelet transform; the ones between the coarsest and the finest are all curvelet transforms. This option only affects the finest scale.

dtypeDTypeLike, optional

dtype of the transform.

class curvelops.curvelops.FDCT3D(dims, axes=(-3, -2, -1), nbscales=None, nbangles_coarse=16, allcurvelets=True, dtype='complex128')[source]#

Bases: FDCT

3D dimensional Curvelet operator. Apply 3D Curvelet Transform along two axes of a multi-dimensional array of size dims.

Parameters#

dimstuple

Number of samples for each dimension.

axestuple, optional

Axes along which FDCT is applied.

nbscalesint, optional

Number of scales (including the coarsest level); Defaults to ceil(log2(min(input_dims)) - 3).

nbangles_coarseint, optional

Number of angles at 2nd coarsest scale.

allcurveletsbool, optional

Use curvelets at the finest (last) scale. If False, a wavelet transform will be used for the finest scale. The coarsest scale is always a wavelet transform; the ones between the coarsest and the finest are all curvelet transforms. This option only affects the finest scale.

dtypeDTypeLike, optional

dtype of the transform.

plot#

curvelops.plot#

Auxiliary functions for plotting.

curvelops.plot.create_axes_grid(rows, cols, kwargs_figure=None, kwargs_gridspec=None, kwargs_subplots=None)[source]#

Creates a grid of axes.

Parameters#

rowsint

Number of rows.

colsint

Number of columns.

kwargs_figureOptional[dict], optional

Arguments to be passed to matplotlib.pyplot.figure.

kwargs_gridspecOptional[dict], optional

Arguments to be passed to matplotlib.gridspec.GridSpec.

kwargs_subplotsOptional[dict], optional

Arguments to be passed to matplotlib.figure.Figure.add_subplot.

Returns#

Tuple[Figure, NDArray]

fig : Figure.

axs : Array of Axes shaped (rows, cols).

Examples#

>>> from curvelops.plot import create_axes_grid
>>> rows, cols = 2, 3
>>> fig, axs = create_axes_grid(
>>>     rows,
>>>     cols,
>>>     kwargs_figure=dict(figsize=(8, 8)),
>>>     kwargs_gridspec=dict(wspace=0.3, hspace=0.3),
>>> )
>>> for irow in range(rows):
>>>     for icol in range(cols):
>>>         axs[irow][icol].plot(np.cos((2 + irow + icol**2) * np.linspace(0, 1)))
>>>         axs[irow][icol].set(title=f"Row, Col: ({irow}, {icol})")
curvelops.plot.create_colorbar(im, ax, size=0.05, pad=0.1, orientation='vertical')[source]#

Create a colorbar.

Divides axis and attaches a colorbar to it.

Parameters#

imAxesImage

Image from which the colorbar will be created. Commonly the output of matplotlib.pyplot.imshow.

axAxes

Axis which to split.

sizefloat, optional

Size of split, by default 0.05. Effectively sets the size of the colorbar.

padfloat, optional`

Padding between colorbar axis and input axis, by default 0.1.

orientationstr, optional

Orientation of the colorbar, by default “vertical”.

Returns#

Tuple[Axes, Colorbar]

cax : Colorbar axis.

cb : Colorbar.

Examples#

>>> import matplotlib.pyplot as plt
>>> from matplotlib.ticker import MultipleLocator
>>> from curvelops.plot import create_colorbar
>>> fig, ax = plt.subplots()
>>> im = ax.imshow([[0]], vmin=-1, vmax=1, cmap="gray")
>>> cax, cb = create_colorbar(im, ax)
>>> cax.yaxis.set_major_locator(MultipleLocator(0.1))
>>> print(cb.vmin)
-1.0
curvelops.plot.create_inset_axes_grid(ax, rows, cols, height=0.5, width=0.5, kwargs_inset_axes=None)[source]#

Create a grid of insets.

The input axis will be overlaid with a grid of insets. Numbering of the axes is top to bottom (rows) and left to right (cols).

Parameters#

axAxes

Input axis.

rowsint

Number of rows.

colsint

Number of columns.

widthfloat, optional

Width of each axis, as a percentage of cols, by default 0.5.

heightfloat, optional

Height of each axis, as a percentage of rows, by default 0.5.

kwargs_inset_axesOptional[dict], optional

Arguments to be passed to matplotlib.axes.Axes.inset_axes.

Returns#

NDArray

Array of Axes shaped (rows, cols).

Examples#

>>> import matplotlib.pyplot as plt
>>> import numpy as np
>>> from curvelops.plot import create_inset_axes_grid
>>> fig, ax = plt.subplots(figsize=(6, 6))
>>> ax.imshow([[0]], extent=[-2, 2, 2, -2], vmin=-1, vmax=1, cmap="gray")
>>> rows, cols = 2, 3
>>> inset_axes = create_inset_axes_grid(
>>>     ax,
>>>     rows,
>>>     cols,
>>>     width=0.5,
>>>     height=0.5,
>>>     kwargs_inset_axes=dict(projection="polar"),
>>> )
>>> nscales = 4
>>> lw = 0.1
>>> for irow in range(rows):
>>>     for icol in range(cols):
>>>         for iscale in range(1, nscales):
>>>             inset_axes[irow][icol].bar(
>>>                 x=0,
>>>                 height=lw,
>>>                 width=2 * np.pi,
>>>                 bottom=((iscale + 1) - 0.5 * lw) / (nscales - 1),
>>>                 color="r",
>>>             )
>>>             inset_axes[irow][icol].set(title=f"Row, Col: ({irow}, {icol})")
>>>             inset_axes[irow][icol].axis("off")
curvelops.plot.curveshow(c_struct, k_space=False, basesize=5, showaxis=False, real=True, kwargs_imshow=None)[source]#

Display curvelet coefficients in each wedge as images.

For each curvelet scale, display a figure with each wedge plotted as an image in its own axis.

Parameters#

c_structFDCTStructLike

Curvelet structure.

k_spacebool, optional

Show curvelet coefficient (False) or its 2D FFT transform (True), by default False.

basesizeint, optional

Base fize of figure, by default 5. Each figure will be sized (basesize * cols, basesize * rows), where rows = floor(sqrt(nangles)) and cols = ceil(nangles / rows)

showaxisbool, optional

Turn on axis lines and labels, by default False.

realbool, optional

Plot real or imaginary part of curvelet coefficients. Only applicable when k_space is False.

kwargs_imshowOptional[dict], optional

Arguments to be passed to matplotlib.pyplot.imshow.

Examples#

>>> import numpy as np
>>> from curvelops import FDCT2D
>>> from curvelops.utils import apply_along_wedges, energy
>>> from curvelops.plot import curveshow
>>> d = np.random.randn(101, 101)
>>> C = FDCT2D(d.shape, nbscales=2, nbangles_coarse=8)
>>> y = C.struct(C @ d)
>>> y_norm = apply_along_wedges(y, lambda w, *_: w / energy(w))
>>> curveshow(
>>>     y_norm,
>>>     basesize=2,
>>>     kwargs_imshow=dict(aspect="auto", vmin=-1, vmax=1, cmap="RdBu")
>>> )

Returns#

List[Figure]

One figure per scale.

curvelops.plot.overlay_arrows(vectors, ax, arrowprops=None)[source]#

Overlay arrows on an axis.

Parameters#

vectorsNDArray

Array shaped (rows, cols, 2), corresponding to a 2D vector field.

axAxes

Axis on which to overlay the arrows.

arrowpropsOptional[dict], optional

Arrow properties, to be passed to matplotlib.pyplot.annotate. By default will be set to dict(facecolor="black", shrink=0.05).

Examples#

>>> import matplotlib.pyplot as plt
>>> import numpy as np
>>> from curvelops.plot import overlay_arrows
>>> fig, ax = plt.subplots(figsize=(8, 10))
>>> ax.imshow([[0]], vmin=-1, vmax=1, extent=[0, 1, 1, 0], cmap="gray")
>>> rows, cols = 3, 4
>>> kvecs = np.array(
>>>     [
>>>         [(1 + x, x * y) for x in (0.5 + np.arange(cols)) / cols]
>>>         for y in (0.5 + np.arange(rows)) / rows
>>>     ]
>>> )
>>> overlay_arrows(
>>>     0.05 * kvecs,
>>>     ax,
>>>     arrowprops=dict(
>>>         facecolor="r",
>>>         shrink=0.05,
>>>         width=10 / cols,
>>>         headwidth=10,
>>>         headlength=10,
>>>     ),
>>> )
curvelops.plot.overlay_disks(c_struct, axes, linewidth=0.5, linecolor='r', map_cmap=True, cmap='gray_r', alpha=1.0, pclip=1.0, map_alpha=False, min_alpha=0.05, normalize='all', annotate=False)[source]#

Overlay curvelet disks over a 2D grid of axes.

Its intended usage is to display the strength of curvelet coefficients of a certain image with a disk display. Given an axes 2D array, each curvelet wedge will be split into rows, cols = axes.shape sub-wedges. The energy of each of these sub-wedges will be mapped to a colormap color and/or transparency.

See Also#

energy_split: Splits a wedge into (rows, cols) wedges and computes the energy of each of these subdivisions.

create_inset_axes_grid: Create a grid of insets.

create_axes_grid: Creates a grid of axes.

curveshow: Display curvelet coefficients in each wedge as images.

Parameters#

c_structFDCTStructLike:

Curvelet coefficients of underlying image.

axesNDArray

2D grid of axes for which disks will be computed.

linewidthfloat, optional

Width of line separating scales, by default 0.5. Will be scaled by 0.1 / nscales internally. Set to zero to disable.

linecolorstr, optional

Color of line separating scales, by default “r”.

map_cmapbool, optional

When enabled, energy will be mapped to the colormap, by default True.

cmapUnion[str, Colormap], optional

Colormap, by default "gray_r".

alphafloat, optional

When using map_cmap, sets a transparecy for all wedges. Has no effect when map_alpha is enabled. By default 1.0.

pclipfloat, optional

Clips the maximum amplitude by this percentage. By default 1.0. Should be between 0.0 and 1.0.

map_alphabool, optional

When enabled, energy will be mapped to the transparency, by default False.

min_alphafloat, optional

When using map_alpha, sets a minimum transparency value. Has no effect when map_alpha is disabled. By default 0.05.

normalizestr, optional

Normalize wedges by:

  • "all" (default)

    Colormap/alpha value of 1.0 will correspond to the maximum energy found across all wedges

  • "scale"

    Colormap/alpha value of 1.0 will correspond to the maximum energy found across all wedges in the same scale.

annotatebool, optional

When true, will display in the middle of the wedge a pair of numbers iscale, iwedge, the index of that scale and that wedge, both starting from zero. This option is useful to understand which directions each wedge corresponds to. By default False.

Examples#

>>> import matplotlib.pyplot as plt
>>> import numpy as np
>>> from curvelops import FDCT2D
>>> from curvelops.utils import apply_along_wedges
>>> from curvelops.plot import create_axes_grid, overlay_disks
>>> x = np.random.randn(50, 100)
>>> C = FDCT2D(x.shape, nbscales=4, nbangles_coarse=8)
>>> y = C.struct(C @ x)
>>> y_ones = apply_along_wedges(y, lambda w, *_: np.ones_like(w))
>>> fig, axes = create_axes_grid(
>>>     1,
>>>     1,
>>>     kwargs_subplots=dict(projection="polar"),
>>>     kwargs_figure=dict(figsize=(8, 8)),
>>> )
>>> overlay_disks(y_ones, axes, annotate=True, cmap="gray")
>>> import matplotlib as mpl
>>> import matplotlib.pyplot as plt
>>> import numpy as np
>>> from mpl_toolkits.axes_grid1 import make_axes_locatable
>>> from curvelops import FDCT2D
>>> from curvelops.plot import create_inset_axes_grid, overlay_disks
>>> from curvelops.utils import apply_along_wedges
>>> plt.rcParams.update({"image.interpolation": "blackman"})
>>> # Construct signal
>>> xlim = [-1.0, 1.0]
>>> ylim = [-0.5, 0.5]
>>> x = np.linspace(*xlim, 201)
>>> z = np.linspace(*ylim, 101)
>>> xm, zm = np.meshgrid(x, z, indexing="ij")
>>> freq = 5
>>> d = np.cos(2 * np.pi * freq * (xm + np.cos(xm) * zm) ** 3)
>>> # Compute curvelet coefficients
>>> C = FDCT2D(d.shape, nbangles_coarse=8, allcurvelets=False)
>>> d_c = C.struct(C @ d)
>>> # Plot original signal
>>> fig, ax = plt.subplots(figsize=(8, 4    ))
>>> ax.imshow(d.T, extent=[*xlim, *(ylim[::-1])], cmap="RdYlBu", vmin=-1, vmax=1)
>>> ax.axis("off")
>>> # Overlay disks
>>> rows, cols = 3, 6
>>> axesin = create_inset_axes_grid(
>>>     ax, rows, cols, width=0.75, kwargs_inset_axes=dict(projection="polar")
>>> )
>>> pclip = 0.2
>>> cmap = "gray_r"
>>> overlay_disks(d_c, axesin, linewidth=0.0, pclip=pclip, cmap=cmap)
>>> # Display disk colorbar
>>> divider = make_axes_locatable(ax)
>>> cax = divider.append_axes("right", size="5%", pad=0.1)
>>> mpl.colorbar.ColorbarBase(
>>>     cax, cmap=cmap, norm=mpl.colors.Normalize(vmin=0, vmax=pclip)
>>> )

typing#

curvelops.typing#

Typing submodule.

utils#

curvelops.utils#

Utility functions for processing curvelets.

curvelops.utils.apply_along_wedges(c_struct, fun)[source]#

Applies a function to each individual wedge.

Parameters#

c_structFDCTStructLike

Input curvelet coefficients in struct format.

funCallable[[NDArray, int, int, int, int], T]

Function to apply to each individual wedge. The function’s arguments are respectively: wedge, wedge index in scale, scale index, number of wedges in scale, number of scales.

Returns#

List[List[T]]

Struct containing the result of applying fun to each wedge.

Examples#

>>> import numpy as np
>>> from curvelops import FDCT2D
>>> from curvelops.utils import apply_along_wedges
>>> x = np.zeros((32, 32))
>>> C = FDCT2D(x.shape, nbscales=3, nbangles_coarse=8, allcurvelets=False)
>>> y = C.struct(C @ x)
>>> apply_along_wedges(y, lambda w, *_: w.shape)
[[(11, 11)],
 [(23, 11),
  (23, 11),
  (11, 23),
  (11, 23),
  (23, 11),
  (23, 11),
  (11, 23),
  (11, 23)],
 [(32, 32)]]
curvelops.utils.array_split_nd(ary, *args)[source]#

Split an array into multiple sub-arrays recursively, possibly unevenly.

See Also#

numpy.array_split : Split an array into multiple sub-arrays.

split_nd: Evenly split an array into multiple sub-arrays recursively.

Parameters#

aryNDArray

Input array.

argsint, optional

Number of splits for each axis of ary. Axis 0 will be split into args[0] subarrays, axis 1 will be into args[1] subarrays, etc. An axis of length l = ary.shape[axis] that should be split into n = args[axis] sections, will return l % n sub-arrays of size l//n + 1 and the rest of size l//n.

Returns#

RecursiveListNDArray

Recursive lists of lists of NDArray. The number of recursions is equivalent to the number arguments in args.

Examples#

>>> from curvelops.utils import array_split_nd
>>> ary = np.outer(1 + np.arange(2), 2 + np.arange(3))
array([[2, 3, 4],
       [4, 6, 8]])
>>> array_split_nd(ary, 2, 3)
[[array([[2]]), array([[3]]), array([[4]])],
 [array([[4]]), array([[6]]), array([[8]])]]
>>> from curvelops.utils import array_split_nd
>>> ary = np.outer(np.arange(3), np.arange(5))
>>> array_split_nd(ary, 2, 3)
[[array([[0, 0],
         [0, 1]]),
  array([[0, 0],
         [2, 3]]),
  array([[0],
         [4]])],
 [array([[0, 2]]), array([[4, 6]]), array([[8]])]]
curvelops.utils.energy(ary)[source]#

Computes the energy of an n-dimensional wedge.

The energy of a vector (flattened n-dimensional array) \((a_0,\ldots,a_{N-1})\) is defined as

\[\sqrt{\frac{1}{N}\sum\limits_{i=0}^{N-1} |a_i|^2}.\]

Parameters#

aryNDArray

Input wedge.

Returns#

float

Energy.

curvelops.utils.energy_split(ary, rows, cols)[source]#

Splits a wedge into (rows, cols) wedges and computes the energy of each of these subdivisions.

See Also#

energy : Computes the energy of a wedge.

Parameters#

aryNDArray

Input wedge.

rowsint

Split axis 0 into rows subdivisions.

colsint

Split axis 1 into cols subdivisions.

Returns#

NDArray

Matrix of shape (rows, cols) containing the energy of each subdivision of the input wedge.

curvelops.utils.ndargmax(ary)[source]#

N-dimensional argmax of array.

Parameters#

aryNDArray

Input array

Examples#

>>> import numpy as np
>>> from curvelops.utils import ndargmax
>>> x = np.zeros((10, 10, 10))
>>> x[1, 1, 1] = 1.0
>>> ndargmax(x)
(1, 1, 1)

Returns#

tuple

N-dimensional index of the maximum of ary.

curvelops.utils.split_nd(ary, *args)[source]#

Evenly split an array into multiple sub-arrays recursively.

See Also#

numpy.split : Split an array into multiple sub-arrays.

array_split_nd: Split an array into multiple sub-arrays recursively, possibly unevenly.

Parameters#

aryNDArray

Input array.

argsint, optional

Number of splits for each axis of ary. Axis 0 will be split into args[0] subarrays, axis 1 will be into args[1] subarrays, etc. If the split cannot be made even for all dimensions, raises an error.

Returns#

RecursiveListNDArray

Recursive lists of lists of NDArray. The number of recursions is equivalent to the number arguments in args.

Examples#

>>> from curvelops.utils import split_nd
>>> ary = np.outer(1 + np.arange(2), 2 + np.arange(3))
array([[2, 3, 4],
       [4, 6, 8]])
>>> split_nd(ary, 2, 3)
[[array([[2]]), array([[3]]), array([[4]])],
 [array([[4]]), array([[6]]), array([[8]])]]
>>> from curvelops.utils import split_nd
>>> ary = np.outer(np.arange(3), np.arange(5))
>>> split_nd(ary, 2, 3)
ValueError: array split does not result in an equal division