API reference

This page provides an auto-generated summary of xrft’s API. For more details and examples, refer to the relevant chapters in the main part of the documentation.

Note

None of xrft’s functions will work correctly in the presence of NaNs or missing data. It’s the user’s responsibility to ensure data are free of NaN or that NaNs have been filled somehow.

xrft

xrft.xrft.cross_phase(da1, da2, dim=None, true_phase=True, **kwargs)[source]

Calculates the cross-phase between da1 and da2.

Returned values are in [-pi, pi].

\[da1' = da1 - \overline{da1};\ \ da2' = da2 - \overline{da2}\]
\[cp = ext{Arg} [\mathbb{F}(da1')^*, \mathbb{F}(da2')]\]
Parameters
da1xarray.DataArray

The data to be transformed

da2xarray.DataArray

The data to be transformed

dimstr or sequence of str, optional

The dimensions along which to take the transformation. If None, all dimensions will be transformed.

true_phaseboolean

If True, the phase information is retained. Set explicitly true_phase = False in cross_spectrum arguments list to ensure future compatibility with numpy-like behavior where the coordinates are disregarded.

kwargsdictsee xrft.fft for argument list
xrft.xrft.cross_spectrum(da1, da2, dim=None, real_dim=None, scaling='density', window_correction=False, true_phase=True, **kwargs)[source]

Calculates the cross spectra of da1 and da2.

\[da1' = da1 - \overline{da1};\ \ da2' = da2 - \overline{da2}\]
\[cs = \mathbb{F}(da1') {\mathbb{F}(da2')}^*\]
Parameters
da1xarray.DataArray

The data to be transformed

da2xarray.DataArray

The data to be transformed

dimstr or sequence of str, optional

The dimensions along which to take the transformation. If None, all dimensions will be transformed.

real_dimstr, optional

Real Fourier transform will be taken along this dimension.

scalingstr, optional

If ‘density’, it will normalize the output to power spectral density If ‘spectrum’, it will normalize the output to power spectrum

window_correctionboolean

If True, it will correct for the energy reduction resulting from applying a non-uniform window. This is the default behaviour of many tools for computing power spectrum (e.g scipy.signal.welch and scipy.signal.periodogram). If scaling = ‘spectrum’, correct the amplitude of peaks in the spectrum. This ensures, for example, that the peak in the one-sided power spectrum of a 10 Hz sine wave with RMS**2 = 10 has a magnitude of 10. If scaling = ‘density’, correct for the energy (integral) of the spectrum. This ensures, for example, that the power spectral density integrates to the square of the RMS of the signal (ie that Parseval’s theorem is satisfied). Note that in most cases, Parseval’s theorem will only be approximately satisfied with this correction as it assumes that the signal being windowed is independent of the window. The correction becomes more accurate as the width of the window gets large in comparison with any noticeable period in the signal. If False, the spectrum gives a representation of the power in the windowed signal. Note that when True, Parseval’s theorem may only be approximately satisfied.

true_phaseboolean

If True, the phase information is retained. Set explicitly true_phase = False in cross_spectrum arguments list to ensure future compatibility with numpy-like behavior where the coordinates are disregarded.

kwargsdictsee xrft.fft for argument list
xrft.xrft.dft(da, dim=None, true_phase=False, true_amplitude=False, **kwargs)[source]

Deprecated function. See fft doc

xrft.xrft.fft(da, spacing_tol=0.001, dim=None, real_dim=None, shift=True, detrend=None, window=None, true_phase=True, true_amplitude=True, chunks_to_segments=False, prefix='freq_', **kwargs)[source]

Perform discrete Fourier transform of xarray data-array da along the specified dimensions.

\[daft = \mathbb{F}(da - \overline{da})\]
Parameters
daxarray.DataArray

The data to be transformed

spacing_tol: float, optional

Spacing tolerance. Fourier transform should not be applied to uneven grid but this restriction can be relaxed with this setting. Use caution.

dimstr or sequence of str, optional

The dimensions along which to take the transformation. If None, all dimensions will be transformed. If the inputs are dask arrays, the arrays must not be chunked along these dimensions.

real_dimstr, optional

Real Fourier transform will be taken along this dimension.

shiftbool, default

Whether to shift the fft output. Default is True, unless real_dim is not None, in which case shift will be set to False always.

detrend{None, ‘constant’, ‘linear’}

If constant, the mean across the transform dimensions will be subtracted before calculating the Fourier transform (FT). If linear, the linear least-square fit will be subtracted before the FT. For linear, only dims of length 1 and 2 are supported.

windowstr, optional

Whether to apply a window to the data before the Fourier transform is taken. A window will be applied to all the dimensions in dim. Please follow scipy.signal.windows’ naming convention.

true_phasebool, optional

If set to False, standard fft algorithm is applied on signal without consideration of coordinates. If set to True, coordinates location are correctly taken into account to evaluate Fourier Tranforrm phase and fftshift is applied on input signal prior to fft (fft algorithm intrinsically considers that input signal is on fftshifted grid).

true_amplitudebool, optional

If set to True, output is multiplied by the spacing of the transformed variables to match theoretical FT amplitude. If set to False, amplitude regularisation by spacing is not applied (as in numpy.fft)

chunks_to_segmentsbool, optional

Whether the data is chunked along the axis to take FFT.

prefixstr

The prefix for the new transformed dimensions.

Returns
daftxarray.DataArray

The output of the Fourier transformation, with appropriate dimensions.

xrft.xrft.fit_loglog(x, y)[source]

Fit a line to isotropic spectra in log-log space

Parameters
xnumpy.array

Coordinate of the data

ynumpy.array

data

Returns
y_fitnumpy.array

The linear fit

afloat64

Slope of the fit

bfloat64

Intercept of the fit

xrft.xrft.idft(daft, dim=None, true_phase=False, true_amplitude=False, **kwargs)[source]

Deprecated function. See ifft doc

xrft.xrft.ifft(daft, spacing_tol=0.001, dim=None, real_dim=None, shift=True, true_phase=True, true_amplitude=True, chunks_to_segments=False, prefix='freq_', lag=None, **kwargs)[source]

Perform inverse discrete Fourier transform of xarray data-array daft along the specified dimensions.

\[da = \mathbb{F}(daft - \overline{daft})\]
Parameters
daftxarray.DataArray

The data to be transformed

spacing_tol: float, optional

Spacing tolerance. Fourier transform should not be applied to uneven grid but this restriction can be relaxed with this setting. Use caution.

dimstr or sequence of str, optional

The dimensions along which to take the transformation. If None, all dimensions will be transformed.

real_dimstr, optional

Real Fourier transform will be taken along this dimension.

shiftbool, default

Whether to shift the fft output. Default is True.

chunks_to_segmentsbool, optional

Whether the data is chunked along the axis to take FFT.

prefixstr

The prefix for the new transformed dimensions.

true_phasebool, optional

If set to False, standard ifft algorithm is applied on signal without consideration of coordinates order. If set to True, coordinates are correctly taken into account to evaluate Inverse Fourier Tranforrm phase and fftshift is applied on input signal prior to ifft (ifft algorithm intrinsically considers that input signal is on fftshifted grid).

true_amplitudebool, optional

If set to True, output is divided by the spacing of the transformed variables to match theoretical IFT amplitude. If set to False, amplitude regularisation by spacing is not applied (as in numpy.ifft)

lagNone, float or sequence of float and/or None, optional

Output coordinates of transformed dimensions will be shifted by corresponding lag values and correct signal phasing will be preserved if true_phase is set to True. If lag is None (default), ‘direct_lag’ attributes of each dimension is used (or set to zero if not found). If defined, lag must have same length as dim. If lag is a sequence, a None element means that ‘direct_lag’ attribute will be used for the corresponding dimension Manually set lag to zero to get output coordinates centered on zero.

Returns
daxarray.DataArray

The output of the Inverse Fourier transformation, with appropriate dimensions.

xrft.xrft.isotropic_cross_spectrum(da1, da2, spacing_tol=0.001, dim=None, shift=True, detrend=None, scaling='density', window=None, window_correction=False, nfactor=4, truncate=False, **kwargs)[source]

Calculates the isotropic spectrum from the two-dimensional power spectrum by taking the azimuthal average.

\[ext{iso}_{cs} = k_r N^{-1} \sum_{N} (\mathbb{F}(da1') {\mathbb{F}(da2')}^*)\]

where \(N\) is the number of azimuthal bins.

Note: the method is not lazy does trigger computations.

Parameters
da1xarray.DataArray

The data to be transformed

da2xarray.DataArray

The data to be transformed

spacing_tol: float (default)

Spacing tolerance. Fourier transform should not be applied to uneven grid but this restriction can be relaxed with this setting. Use caution.

dimlist (optional)

The dimensions along which to take the transformation. If None, all dimensions will be transformed.

shiftbool (optional)

Whether to shift the fft output.

detrendstr (optional)

If constant, the mean across the transform dimensions will be subtracted before calculating the Fourier transform (FT). If linear, the linear least-square fit will be subtracted before the FT.

densitylist (optional)

If true, it will normalize the spectrum to spectral density

windowstr (optional)

Whether to apply a window to the data before the Fourier transform is taken. Please adhere to scipy.signal.windows for naming convention.

nfactorint (optional)

Ratio of number of bins to take the azimuthal averaging with the data size. Default is 4.

truncatebool, optional

If True, the spectrum will be truncated for wavenumbers larger than the Nyquist wavenumber.

Returns
iso_csxarray.DataArray

Isotropic cross spectrum

xrft.xrft.isotropic_power_spectrum(da, spacing_tol=0.001, dim=None, shift=True, detrend=None, scaling='density', window=None, window_correction=False, nfactor=4, truncate=False, **kwargs)[source]

Calculates the isotropic spectrum from the two-dimensional power spectrum by taking the azimuthal average.

\[ext{iso}_{ps} = k_r N^{-1} \sum_{N} |\mathbb{F}(da')|^2\]

where \(N\) is the number of azimuthal bins.

Note: the method is not lazy does trigger computations.

Parameters
daxarray.DataArray

The data to be transformed

spacing_tol: float, optional

Spacing tolerance. Fourier transform should not be applied to uneven grid but this restriction can be relaxed with this setting. Use caution.

dimlist, optional

The dimensions along which to take the transformation. If None, all dimensions will be transformed.

shiftbool, optional

Whether to shift the fft output.

detrendstr, optional

If constant, the mean across the transform dimensions will be subtracted before calculating the Fourier transform (FT). If linear, the linear least-square fit will be subtracted before the FT.

densitylist, optional

If true, it will normalize the spectrum to spectral density

windowstr, optional

Whether to apply a window to the data before the Fourier transform is taken. Please adhere to scipy.signal.windows for naming convention.

nfactorint, optional

Ratio of number of bins to take the azimuthal averaging with the data size. Default is 4.

truncatebool, optional

If True, the spectrum will be truncated for wavenumbers larger than the Nyquist wavenumber.

Returns
iso_psxarray.DataArray

Isotropic power spectrum

xrft.xrft.isotropize(ps, fftdim, nfactor=4, truncate=True, complx=False)[source]

Isotropize a 2D power spectrum or cross spectrum by taking an azimuthal average.

\[ext{iso}_{ps} = k_r N^{-1} \sum_{N} |\mathbb{F}(da')|^2\]

where \(N\) is the number of azimuthal bins.

Parameters
psxarray.DataArray

The power spectrum or cross spectrum to be isotropized.

fftdimlist

The fft dimensions overwhich the isotropization must be performed.

nfactorint, optional

Ratio of number of bins to take the azimuthal averaging with the data size. Default is 4.

truncatebool, optional

If True, the spectrum will be truncated for wavenumbers larger than the Nyquist wavenumber.

complxbool, optional

If True, isotropize allows for complex numbers.

xrft.xrft.power_spectrum(da, dim=None, real_dim=None, scaling='density', window_correction=False, **kwargs)[source]

Calculates the power spectrum of da.

\[\]

da’ = da - overline{da} .. math:: ps = mathbb{F}(da’) {mathbb{F}(da’)}^*

Parameters
daxarray.DataArray

The data to be transformed

dimstr or sequence of str, optional

The dimensions along which to take the transformation. If None, all dimensions will be transformed.

real_dimstr, optional

Real Fourier transform will be taken along this dimension.

scalingstr, optional

If ‘density’, it will normalize the output to power spectral density If ‘spectrum’, it will normalize the output to power spectrum

window_correctionboolean

If True, it will correct for the energy reduction resulting from applying a non-uniform window. This is the default behaviour of many tools for computing power spectrum (e.g scipy.signal.welch and scipy.signal.periodogram). If scaling = ‘spectrum’, correct the amplitude of peaks in the spectrum. This ensures, for example, that the peak in the one-sided power spectrum of a 10 Hz sine wave with RMS**2 = 10 has a magnitude of 10. If scaling = ‘density’, correct for the energy (integral) of the spectrum. This ensures, for example, that the power spectral density integrates to the square of the RMS of the signal (ie that Parseval’s theorem is satisfied). Note that in most cases, Parseval’s theorem will only be approximately satisfied with this correction as it assumes that the signal being windowed is independent of the window. The correction becomes more accurate as the width of the window gets large in comparison with any noticeable period in the signal. If False, the spectrum gives a representation of the power in the windowed signal. Note that when True, Parseval’s theorem may only be approximately satisfied.

kwargsdictsee xrft.fft for argument list

detrend

You also may wish to use xrft’s detrend function on its own.

Functions for detrending xarray data.

xrft.detrend.detrend(da, dim, detrend_type='constant')[source]

Detrend a DataArray

Parameters
daxarray.DataArray

The data to detrend

dimstr or list

Dimensions along which to apply detrend. Can be either one dimension or a list with two dimensions. Higher-dimensional detrending is not supported. If dask data are passed, the data must be chunked along dim.

detrend_type{‘constant’, ‘linear’}

If constant, a constant offset will be removed from each dim. If linear, a linear least-squares fit will be estimated and removed from the data.

Returns
daxarray.DataArray

The detrended data.

Notes

This function will act lazily in the presence of dask arrays on the input.

padding

Pad and unpad arrays and its coordinates so they can be used for computing FFTs.

Functions to pad and unpad a N-dimensional regular grid

xrft.padding.pad(da, pad_width=None, mode='constant', stat_length=None, constant_values=0, end_values=None, reflect_type=None, **pad_width_kwargs)[source]

Pad array with evenly spaced coordinates

Wraps the xarray.DataArray.pad() method but also pads the evenly spaced coordinates by extrapolation using the same coordinate spacing. The pad_width used for each coordinate is stored as one of its attributes.

Parameters
daxarray.DataArray

Array to be padded. The coordinates along which the array will be padded must be evenly spaced.

pad_widthmapping of hashable to tuple of int

Mapping with the form of {dim: (pad_before, pad_after)} describing the number of values padded along each dimension. {dim: pad} is a shortcut for pad_before = pad_after = pad

modestr, default: “constant”

One of the following string values (taken from numpy docs). - constant: Pads with a constant value. - edge: Pads with the edge values of array. - linear_ramp: Pads with the linear ramp between end_value and the

array edge value.

  • maximum: Pads with the maximum value of all or part of the vector along each axis.

  • mean: Pads with the mean value of all or part of the vector along each axis.

  • median: Pads with the median value of all or part of the vector along each axis.

  • minimum: Pads with the minimum value of all or part of the vector along each axis.

  • reflect: Pads with the reflection of the vector mirrored on the first and last values of the vector along each axis.

  • symmetric: Pads with the reflection of the vector mirrored along the edge of the array.

  • wrap: Pads with the wrap of the vector along the axis. The first values are used to pad the end and the end values are used to pad the beginning.

stat_lengthint, tuple or mapping of hashable to tuple, default: None

Used in ‘maximum’, ‘mean’, ‘median’, and ‘minimum’. Number of values at edge of each axis used to calculate the statistic value. {dim_1: (before_1, after_1), … dim_N: (before_N, after_N)} unique statistic lengths along each dimension. ((before, after),) yields same before and after statistic lengths for each dimension. (stat_length,) or int is a shortcut for before = after = statistic length for all axes. Default is None, to use the entire axis.

constant_valuesscalar, tuple or mapping of hashable to tuple, default: 0

Used in ‘constant’. The values to set the padded values for each axis. {dim_1: (before_1, after_1), ... dim_N: (before_N, after_N)} unique pad constants along each dimension. ((before, after),) yields same before and after constants for each dimension. (constant,) or constant is a shortcut for before = after = constant for all dimensions. Default is 0.

end_valuesscalar, tuple or mapping of hashable to tuple, default: 0

Used in ‘linear_ramp’. The values used for the ending value of the linear_ramp and that will form the edge of the padded array. {dim_1: (before_1, after_1), ... dim_N: (before_N, after_N)} unique end values along each dimension. ((before, after),) yields same before and after end values for each axis. (constant,) or constant is a shortcut for before = after = constant for all axes. Default is 0.

reflect_type{“even”, “odd”}, optional

Used in “reflect”, and “symmetric”. The “even” style is the default with an unaltered reflection around the edge value. For the “odd” style, the extended part of the array is created by subtracting the reflected values from two times the edge value.

**pad_width_kwargs

The keyword arguments form of pad_width. One of pad_width or pad_width_kwargs must be provided.

Returns
da_paddedxarray.DataArray

Examples

>>> import xarray as xr
>>> da = xr.DataArray(
...     [[1, 2, 3], [4, 5, 6], [7, 8, 9]],
...     coords={"x": [0, 1, 2], "y": [-5, -4, -3]},
...     dims=("y", "x"),
... )
>>> da_padded = pad(da, x=2, y=1)
>>> da_padded
<xarray.DataArray (y: 5, x: 7)>
array([[0, 0, 0, 0, 0, 0, 0],
       [0, 0, 1, 2, 3, 0, 0],
       [0, 0, 4, 5, 6, 0, 0],
       [0, 0, 7, 8, 9, 0, 0],
       [0, 0, 0, 0, 0, 0, 0]])
Coordinates:
  * x        (x) int64 -2 -1 0 1 2 3 4
  * y        (y) int64 -6 -5 -4 -3 -2
>>> da_padded.x
<xarray.DataArray 'x' (x: 7)>
array([-2, -1,  0,  1,  2,  3,  4])
Coordinates:
  * x        (x) int64 -2 -1 0 1 2 3 4
Attributes:
    pad_width:  2
>>> da_padded.y
<xarray.DataArray 'y' (y: 5)>
array([-6, -5, -4, -3, -2])
Coordinates:
  * y        (y) int64 -6 -5 -4 -3 -2
Attributes:
    pad_width:  1

Asymmetric padding

>>> da_padded = pad(da, x=(1, 4))
>>> da_padded
<xarray.DataArray (y: 3, x: 8)>
array([[0, 1, 2, 3, 0, 0, 0, 0],
       [0, 4, 5, 6, 0, 0, 0, 0],
       [0, 7, 8, 9, 0, 0, 0, 0]])
Coordinates:
  * x        (x) int64 -1 0 1 2 3 4 5 6
  * y        (y) int64 -5 -4 -3
>>> da_padded.x
<xarray.DataArray 'x' (x: 8)>
array([-1,  0,  1,  2,  3,  4,  5,  6])
Coordinates:
  * x        (x) int64 -1 0 1 2 3 4 5 6
Attributes:
    pad_width:  (1, 4)
xrft.padding.unpad(da, pad_width=None, **pad_width_kwargs)[source]

Unpad an array and its coordinates

Undo the padding process of the xrft.pad() function by slicing the passed xarray.DataArray and its coordinates.

Parameters
daxarray.DataArray

Padded array. The coordinates along which the array will be padded must be evenly spaced.

Returns
da_unpadedxarray.DataArray

Unpadded array.

pad_widthmapping of hashable to tuple of int (optional)

Mapping with the form of {dim: (pad_before, pad_after)} describing the number of values padded along each dimension. {dim: pad} is a shortcut for pad_before = pad_after = pad. If None, then the pad_width for each coordinate is read from their pad_width attribute.

**pad_width_kwargs (optional)

The keyword arguments form of pad_width. Pass pad_width or pad_width_kwargs.

See also

xrft.pad()

Examples

>>> import xarray as xr
>>> da = xr.DataArray(
...     [[1, 2, 3], [4, 5, 6], [7, 8, 9]],
...     coords={"x": [0, 1, 2], "y": [-5, -4, -3]},
...     dims=("y", "x"),
... )
>>> da_padded = pad(da, x=2, y=1)
>>> da_padded
<xarray.DataArray (y: 5, x: 7)>
array([[0, 0, 0, 0, 0, 0, 0],
       [0, 0, 1, 2, 3, 0, 0],
       [0, 0, 4, 5, 6, 0, 0],
       [0, 0, 7, 8, 9, 0, 0],
       [0, 0, 0, 0, 0, 0, 0]])
Coordinates:
  * x        (x) int64 -2 -1 0 1 2 3 4
  * y        (y) int64 -6 -5 -4 -3 -2
>>> unpad(da_padded)
<xarray.DataArray (y: 3, x: 3)>
array([[1, 2, 3],
       [4, 5, 6],
       [7, 8, 9]])
Coordinates:
  * x        (x) int64 0 1 2
  * y        (y) int64 -5 -4 -3

Custom pad_width

>>> unpad(da_padded, x=1, y=1)
<xarray.DataArray (y: 3, x: 5)>
array([[0, 1, 2, 3, 0],
       [0, 4, 5, 6, 0],
       [0, 7, 8, 9, 0]])
Coordinates:
  * x        (x) int64 -1 0 1 2 3
  * y        (y) int64 -5 -4 -3