Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Deprecating dft and replacing it with fft #161

Merged
merged 3 commits into from
Aug 24, 2021
Merged
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
90 changes: 44 additions & 46 deletions xrft/xrft.py
Original file line number Diff line number Diff line change
Expand Up @@ -237,47 +237,39 @@ def _lag_coord(coord):
return lag.data


def fft(da, dim=None, **kwargs):
def dft(
da, dim=None, true_phase=False, true_amplitude=False, **kwargs
): # pragma: no cover
"""
da : `xarray.DataArray`
The data to be transformed
dim : str 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.
kwargs: See xrft.dft for argument list
Deprecated function. See fft doc
"""
if kwargs.pop("true_phase", False):
warnings.warn("true_phase argument is ignored in xrft.fft")
if kwargs.pop("true_amplitude", False):
warnings.warn("true_amplitude argument is ignored in xrft.fft")
with warnings.catch_warnings():
warnings.simplefilter("ignore")
return dft(da, dim=dim, true_phase=False, true_amplitude=False, **kwargs)
msg = (
"This function has been renamed and will disappear in the future."
+ " Please use `fft` instead"
)
warnings.warn(msg, FutureWarning)
return fft(
da, dim=dim, true_phase=true_phase, true_amplitude=true_amplitude, **kwargs
)


def ifft(daft, dim=None, **kwargs):
def idft(
daft, dim=None, true_phase=False, true_amplitude=False, **kwargs
): # pragma: no cover
"""
daft : `xarray.DataArray`
The data to be transformed
dim : str 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.
kwargs: See xrft.idft for argument list
Deprecated function. See ifft doc
"""
if kwargs.pop("true_phase", False):
warnings.warn("true_phase argument is ignored in xrft.ifft")
if kwargs.pop("true_amplitude", False):
warnings.warn("true_amplitude argument is ignored in xrft.ifft")
if kwargs.pop("lag", False):
warnings.warn("lag argument is ignored in xrft.ifft")
with warnings.catch_warnings():
warnings.simplefilter("ignore")
return idft(daft, dim=dim, true_phase=False, true_amplitude=False, **kwargs)
msg = (
"This function has been renamed and will disappear in the future."
+ " Please use `ifft` instead"
)
warnings.warn(msg, FutureWarning)
return ifft(
daft, dim=dim, true_phase=true_phase, true_amplitude=true_amplitude, **kwargs
)


def dft(
def fft(
da,
spacing_tol=1e-3,
dim=None,
Expand Down Expand Up @@ -376,13 +368,13 @@ def dft(
*[d for d in da.dims if d not in [real_dim]] + [real_dim]
) # dimension for real transformed is moved at the end

fft = _fft_module(da)
fftm = _fft_module(da)

if real_dim is None:
fft_fn = fft.fftn
fft_fn = fftm.fftn
else:
shift = False
fft_fn = fft.rfftn
fft_fn = fftm.rfftn

# the axes along which to take ffts
axis_num = [
Expand Down Expand Up @@ -426,13 +418,14 @@ def dft(
da.get_axis_num(d) for d in dim if da[d][-1] < da[d][0]
] # handling decreasing coordinates
f = fft_fn(
fft.ifftshift(np.flip(da, axis=reversed_axis), axes=axis_num), axes=axis_num
fftm.ifftshift(np.flip(da, axis=reversed_axis), axes=axis_num),
axes=axis_num,
)
else:
f = fft_fn(da.data, axes=axis_num)

if shift:
f = fft.fftshift(f, axes=axis_num)
f = fftm.fftshift(f, axes=axis_num)

k = _freq(N, delta_x, real_dim, shift)

Expand Down Expand Up @@ -464,7 +457,7 @@ def dft(
) # Do nothing if da was not transposed


def idft(
def ifft(
daft,
spacing_tol=1e-3,
dim=None,
Expand Down Expand Up @@ -716,7 +709,7 @@ def power_spectrum(
{"true_amplitude": True, "true_phase": False}
) # true_phase do not matter in power_spectrum

daft = dft(da, dim=dim, real_dim=real_dim, **kwargs)
daft = fft(da, dim=dim, real_dim=real_dim, **kwargs)
updated_dims = [
d for d in daft.dims if (d not in da.dims and "segment" not in d)
] # Transformed dimensions
Expand Down Expand Up @@ -769,6 +762,7 @@ def cross_spectrum(
real_dim=None,
scaling="density",
window_correction=False,
true_phase=False,
**kwargs,
):
"""
Expand Down Expand Up @@ -803,7 +797,7 @@ def cross_spectrum(
kwargs : dict : see xrft.dft for argument list
"""

if "true_phase" not in kwargs:
if not true_phase:
msg = (
"true_phase flag will be set to True in future version of xrft.dft possibly impacting cross_spectrum output. "
+ "Set explicitely true_phase = False in cross_spectrum arguments list to ensure future compatibility "
Expand All @@ -829,8 +823,8 @@ def cross_spectrum(

kwargs.update({"true_amplitude": True})

daft1 = dft(da1, dim=dim, real_dim=real_dim, **kwargs)
daft2 = dft(da2, dim=dim, real_dim=real_dim, **kwargs)
daft1 = fft(da1, dim=dim, real_dim=real_dim, true_phase=true_phase, **kwargs)
daft2 = fft(da2, dim=dim, real_dim=real_dim, true_phase=true_phase, **kwargs)

if daft1.dims != daft2.dims:
raise ValueError("The two datasets have different dimensions")
Expand Down Expand Up @@ -880,7 +874,7 @@ def cross_spectrum(
return cs


def cross_phase(da1, da2, dim=None, **kwargs):
def cross_phase(da1, da2, dim=None, true_phase=False, **kwargs):
"""
Calculates the cross-phase between da1 and da2.

Expand All @@ -899,15 +893,17 @@ def cross_phase(da1, da2, dim=None, **kwargs):
The data to be transformed
kwargs : dict : see xrft.dft for argument list
"""
if "true_phase" not in kwargs:
if not true_phase:
msg = (
"true_phase flag will be set to True in future version of xrft.dft possibly impacting cross_phase output. "
+ "Set explicitely true_phase = False in cross_spectrum arguments list to ensure future compatibility "
+ "with numpy-like behavior where the coordinates are disregarded."
)
warnings.warn(msg, FutureWarning)

cp = xr.ufuncs.angle(cross_spectrum(da1, da2, dim=dim, **kwargs))
cp = xr.ufuncs.angle(
cross_spectrum(da1, da2, dim=dim, true_phase=true_phase, **kwargs)
)

if da1.name and da2.name:
cp.name = "{}_{}_phase".format(da1.name, da2.name)
Expand Down Expand Up @@ -1134,6 +1130,7 @@ def isotropic_power_spectrum(
scaling=scaling,
window_correction=window_correction,
window=window,
**kwargs,
)

fftdim = ["freq_" + d for d in dim]
Expand Down Expand Up @@ -1237,6 +1234,7 @@ def isotropic_cross_spectrum(
scaling=scaling,
window_correction=window_correction,
window=window,
**kwargs,
)

fftdim = ["freq_" + d for d in dim]
Expand Down