# Copyright (C) 2003-2005 Peter J. Verveer # # Redistribution and use in source and binary forms, with or without # modification, are permitted provided that the following conditions # are met: # # 1. Redistributions of source code must retain the above copyright # notice, this list of conditions and the following disclaimer. # # 2. Redistributions in binary form must reproduce the above # copyright notice, this list of conditions and the following # disclaimer in the documentation and/or other materials provided # with the distribution. # # 3. The name of the author may not be used to endorse or promote # products derived from this software without specific prior # written permission. # # THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS # OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED # WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE # ARE DISCLAIMED. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY # DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL # DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE # GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS # INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, # WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING # NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS # SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. import itertools import warnings import numpy from scipy._lib._util import normalize_axis_index from scipy import special from . import _ni_support from . import _nd_image from ._ni_docstrings import docfiller __all__ = ['spline_filter1d', 'spline_filter', 'geometric_transform', 'map_coordinates', 'affine_transform', 'shift', 'zoom', 'rotate'] @docfiller def spline_filter1d(input, order=3, axis=-1, output=numpy.float64, mode='mirror'): """ Calculate a 1-D spline filter along the given axis. The lines of the array along the given axis are filtered by a spline filter. The order of the spline must be >= 2 and <= 5. Parameters ---------- %(input)s order : int, optional The order of the spline, default is 3. axis : int, optional The axis along which the spline filter is applied. Default is the last axis. output : ndarray or dtype, optional The array in which to place the output, or the dtype of the returned array. Default is ``numpy.float64``. %(mode_interp_mirror)s Returns ------- spline_filter1d : ndarray The filtered input. See Also -------- spline_filter : Multidimensional spline filter. Notes ----- All of the interpolation functions in `ndimage` do spline interpolation of the input image. If using B-splines of `order > 1`, the input image values have to be converted to B-spline coefficients first, which is done by applying this 1-D filter sequentially along all axes of the input. All functions that require B-spline coefficients will automatically filter their inputs, a behavior controllable with the `prefilter` keyword argument. For functions that accept a `mode` parameter, the result will only be correct if it matches the `mode` used when filtering. For complex-valued `input`, this function processes the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- We can filter an image using 1-D spline along the given axis: >>> from scipy.ndimage import spline_filter1d >>> import numpy as np >>> import matplotlib.pyplot as plt >>> orig_img = np.eye(20) # create an image >>> orig_img[10, :] = 1.0 >>> sp_filter_axis_0 = spline_filter1d(orig_img, axis=0) >>> sp_filter_axis_1 = spline_filter1d(orig_img, axis=1) >>> f, ax = plt.subplots(1, 3, sharex=True) >>> for ind, data in enumerate([[orig_img, "original image"], ... [sp_filter_axis_0, "spline filter (axis=0)"], ... [sp_filter_axis_1, "spline filter (axis=1)"]]): ... ax[ind].imshow(data[0], cmap='gray_r') ... ax[ind].set_title(data[1]) >>> plt.tight_layout() >>> plt.show() """ if order < 0 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, complex_output=complex_output) if complex_output: spline_filter1d(input.real, order, axis, output.real, mode) spline_filter1d(input.imag, order, axis, output.imag, mode) return output if order in [0, 1]: output[...] = numpy.array(input) else: mode = _ni_support._extend_mode_to_code(mode) axis = normalize_axis_index(axis, input.ndim) _nd_image.spline_filter1d(input, order, axis, output, mode) return output @docfiller def spline_filter(input, order=3, output=numpy.float64, mode='mirror'): """ Multidimensional spline filter. Parameters ---------- %(input)s order : int, optional The order of the spline, default is 3. output : ndarray or dtype, optional The array in which to place the output, or the dtype of the returned array. Default is ``numpy.float64``. %(mode_interp_mirror)s Returns ------- spline_filter : ndarray Filtered array. Has the same shape as `input`. See Also -------- spline_filter1d : Calculate a 1-D spline filter along the given axis. Notes ----- The multidimensional filter is implemented as a sequence of 1-D spline filters. The intermediate arrays are stored in the same data type as the output. Therefore, for output types with a limited precision, the results may be imprecise because intermediate results may be stored with insufficient precision. For complex-valued `input`, this function processes the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- We can filter an image using multidimentional splines: >>> from scipy.ndimage import spline_filter >>> import numpy as np >>> import matplotlib.pyplot as plt >>> orig_img = np.eye(20) # create an image >>> orig_img[10, :] = 1.0 >>> sp_filter = spline_filter(orig_img, order=3) >>> f, ax = plt.subplots(1, 2, sharex=True) >>> for ind, data in enumerate([[orig_img, "original image"], ... [sp_filter, "spline filter"]]): ... ax[ind].imshow(data[0], cmap='gray_r') ... ax[ind].set_title(data[1]) >>> plt.tight_layout() >>> plt.show() """ if order < 2 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, complex_output=complex_output) if complex_output: spline_filter(input.real, order, output.real, mode) spline_filter(input.imag, order, output.imag, mode) return output if order not in [0, 1] and input.ndim > 0: for axis in range(input.ndim): spline_filter1d(input, order, axis, output=output, mode=mode) input = output else: output[...] = input[...] return output def _prepad_for_spline_filter(input, mode, cval): if mode in ['nearest', 'grid-constant']: npad = 12 if mode == 'grid-constant': padded = numpy.pad(input, npad, mode='constant', constant_values=cval) elif mode == 'nearest': padded = numpy.pad(input, npad, mode='edge') else: # other modes have exact boundary conditions implemented so # no prepadding is needed npad = 0 padded = input return padded, npad @docfiller def geometric_transform(input, mapping, output_shape=None, output=None, order=3, mode='constant', cval=0.0, prefilter=True, extra_arguments=(), extra_keywords={}): """ Apply an arbitrary geometric transform. The given mapping function is used to find, for each point in the output, the corresponding coordinates in the input. The value of the input at those coordinates is determined by spline interpolation of the requested order. Parameters ---------- %(input)s mapping : {callable, scipy.LowLevelCallable} A callable object that accepts a tuple of length equal to the output array rank, and returns the corresponding input coordinates as a tuple of length equal to the input array rank. output_shape : tuple of ints, optional Shape tuple. %(output)s order : int, optional The order of the spline interpolation, default is 3. The order has to be in the range 0-5. %(mode_interp_constant)s %(cval)s %(prefilter)s extra_arguments : tuple, optional Extra arguments passed to `mapping`. extra_keywords : dict, optional Extra keywords passed to `mapping`. Returns ------- output : ndarray The filtered input. See Also -------- map_coordinates, affine_transform, spline_filter1d Notes ----- This function also accepts low-level callback functions with one the following signatures and wrapped in `scipy.LowLevelCallable`: .. code:: c int mapping(npy_intp *output_coordinates, double *input_coordinates, int output_rank, int input_rank, void *user_data) int mapping(intptr_t *output_coordinates, double *input_coordinates, int output_rank, int input_rank, void *user_data) The calling function iterates over the elements of the output array, calling the callback function at each element. The coordinates of the current output element are passed through ``output_coordinates``. The callback function must return the coordinates at which the input must be interpolated in ``input_coordinates``. The rank of the input and output arrays are given by ``input_rank`` and ``output_rank`` respectively. ``user_data`` is the data pointer provided to `scipy.LowLevelCallable` as-is. The callback function must return an integer error status that is zero if something went wrong and one otherwise. If an error occurs, you should normally set the Python error status with an informative message before returning, otherwise a default error message is set by the calling function. In addition, some other low-level function pointer specifications are accepted, but these are for backward compatibility only and should not be used in new code. For complex-valued `input`, this function transforms the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- >>> import numpy as np >>> from scipy.ndimage import geometric_transform >>> a = np.arange(12.).reshape((4, 3)) >>> def shift_func(output_coords): ... return (output_coords[0] - 0.5, output_coords[1] - 0.5) ... >>> geometric_transform(a, shift_func) array([[ 0. , 0. , 0. ], [ 0. , 1.362, 2.738], [ 0. , 4.812, 6.187], [ 0. , 8.263, 9.637]]) >>> b = [1, 2, 3, 4, 5] >>> def shift_func(output_coords): ... return (output_coords[0] - 3,) ... >>> geometric_transform(b, shift_func, mode='constant') array([0, 0, 0, 1, 2]) >>> geometric_transform(b, shift_func, mode='nearest') array([1, 1, 1, 1, 2]) >>> geometric_transform(b, shift_func, mode='reflect') array([3, 2, 1, 1, 2]) >>> geometric_transform(b, shift_func, mode='wrap') array([2, 3, 4, 1, 2]) """ if order < 0 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) if output_shape is None: output_shape = input.shape if input.ndim < 1 or len(output_shape) < 1: raise RuntimeError('input and output rank must be > 0') complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, shape=output_shape, complex_output=complex_output) if complex_output: kwargs = dict(order=order, mode=mode, prefilter=prefilter, output_shape=output_shape, extra_arguments=extra_arguments, extra_keywords=extra_keywords) geometric_transform(input.real, mapping, output=output.real, cval=numpy.real(cval), **kwargs) geometric_transform(input.imag, mapping, output=output.imag, cval=numpy.imag(cval), **kwargs) return output if prefilter and order > 1: padded, npad = _prepad_for_spline_filter(input, mode, cval) filtered = spline_filter(padded, order, output=numpy.float64, mode=mode) else: npad = 0 filtered = input mode = _ni_support._extend_mode_to_code(mode) _nd_image.geometric_transform(filtered, mapping, None, None, None, output, order, mode, cval, npad, extra_arguments, extra_keywords) return output @docfiller def map_coordinates(input, coordinates, output=None, order=3, mode='constant', cval=0.0, prefilter=True): """ Map the input array to new coordinates by interpolation. The array of coordinates is used to find, for each point in the output, the corresponding coordinates in the input. The value of the input at those coordinates is determined by spline interpolation of the requested order. The shape of the output is derived from that of the coordinate array by dropping the first axis. The values of the array along the first axis are the coordinates in the input array at which the output value is found. Parameters ---------- %(input)s coordinates : array_like The coordinates at which `input` is evaluated. %(output)s order : int, optional The order of the spline interpolation, default is 3. The order has to be in the range 0-5. %(mode_interp_constant)s %(cval)s %(prefilter)s Returns ------- map_coordinates : ndarray The result of transforming the input. The shape of the output is derived from that of `coordinates` by dropping the first axis. See Also -------- spline_filter, geometric_transform, scipy.interpolate Notes ----- For complex-valued `input`, this function maps the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- >>> from scipy import ndimage >>> import numpy as np >>> a = np.arange(12.).reshape((4, 3)) >>> a array([[ 0., 1., 2.], [ 3., 4., 5.], [ 6., 7., 8.], [ 9., 10., 11.]]) >>> ndimage.map_coordinates(a, [[0.5, 2], [0.5, 1]], order=1) array([ 2., 7.]) Above, the interpolated value of a[0.5, 0.5] gives output[0], while a[2, 1] is output[1]. >>> inds = np.array([[0.5, 2], [0.5, 4]]) >>> ndimage.map_coordinates(a, inds, order=1, cval=-33.3) array([ 2. , -33.3]) >>> ndimage.map_coordinates(a, inds, order=1, mode='nearest') array([ 2., 8.]) >>> ndimage.map_coordinates(a, inds, order=1, cval=0, output=bool) array([ True, False], dtype=bool) """ if order < 0 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) coordinates = numpy.asarray(coordinates) if numpy.iscomplexobj(coordinates): raise TypeError('Complex type not supported') output_shape = coordinates.shape[1:] if input.ndim < 1 or len(output_shape) < 1: raise RuntimeError('input and output rank must be > 0') if coordinates.shape[0] != input.ndim: raise RuntimeError('invalid shape for coordinate array') complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, shape=output_shape, complex_output=complex_output) if complex_output: kwargs = dict(order=order, mode=mode, prefilter=prefilter) map_coordinates(input.real, coordinates, output=output.real, cval=numpy.real(cval), **kwargs) map_coordinates(input.imag, coordinates, output=output.imag, cval=numpy.imag(cval), **kwargs) return output if prefilter and order > 1: padded, npad = _prepad_for_spline_filter(input, mode, cval) filtered = spline_filter(padded, order, output=numpy.float64, mode=mode) else: npad = 0 filtered = input mode = _ni_support._extend_mode_to_code(mode) _nd_image.geometric_transform(filtered, None, coordinates, None, None, output, order, mode, cval, npad, None, None) return output @docfiller def affine_transform(input, matrix, offset=0.0, output_shape=None, output=None, order=3, mode='constant', cval=0.0, prefilter=True): """ Apply an affine transformation. Given an output image pixel index vector ``o``, the pixel value is determined from the input image at position ``np.dot(matrix, o) + offset``. This does 'pull' (or 'backward') resampling, transforming the output space to the input to locate data. Affine transformations are often described in the 'push' (or 'forward') direction, transforming input to output. If you have a matrix for the 'push' transformation, use its inverse (:func:`numpy.linalg.inv`) in this function. Parameters ---------- %(input)s matrix : ndarray The inverse coordinate transformation matrix, mapping output coordinates to input coordinates. If ``ndim`` is the number of dimensions of ``input``, the given matrix must have one of the following shapes: - ``(ndim, ndim)``: the linear transformation matrix for each output coordinate. - ``(ndim,)``: assume that the 2-D transformation matrix is diagonal, with the diagonal specified by the given value. A more efficient algorithm is then used that exploits the separability of the problem. - ``(ndim + 1, ndim + 1)``: assume that the transformation is specified using homogeneous coordinates [1]_. In this case, any value passed to ``offset`` is ignored. - ``(ndim, ndim + 1)``: as above, but the bottom row of a homogeneous transformation matrix is always ``[0, 0, ..., 1]``, and may be omitted. offset : float or sequence, optional The offset into the array where the transform is applied. If a float, `offset` is the same for each axis. If a sequence, `offset` should contain one value for each axis. output_shape : tuple of ints, optional Shape tuple. %(output)s order : int, optional The order of the spline interpolation, default is 3. The order has to be in the range 0-5. %(mode_interp_constant)s %(cval)s %(prefilter)s Returns ------- affine_transform : ndarray The transformed input. Notes ----- The given matrix and offset are used to find for each point in the output the corresponding coordinates in the input by an affine transformation. The value of the input at those coordinates is determined by spline interpolation of the requested order. Points outside the boundaries of the input are filled according to the given mode. .. versionchanged:: 0.18.0 Previously, the exact interpretation of the affine transformation depended on whether the matrix was supplied as a 1-D or a 2-D array. If a 1-D array was supplied to the matrix parameter, the output pixel value at index ``o`` was determined from the input image at position ``matrix * (o + offset)``. For complex-valued `input`, this function transforms the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. References ---------- .. [1] https://en.wikipedia.org/wiki/Homogeneous_coordinates """ if order < 0 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) if output_shape is None: if isinstance(output, numpy.ndarray): output_shape = output.shape else: output_shape = input.shape if input.ndim < 1 or len(output_shape) < 1: raise RuntimeError('input and output rank must be > 0') complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, shape=output_shape, complex_output=complex_output) if complex_output: kwargs = dict(offset=offset, output_shape=output_shape, order=order, mode=mode, prefilter=prefilter) affine_transform(input.real, matrix, output=output.real, cval=numpy.real(cval), **kwargs) affine_transform(input.imag, matrix, output=output.imag, cval=numpy.imag(cval), **kwargs) return output if prefilter and order > 1: padded, npad = _prepad_for_spline_filter(input, mode, cval) filtered = spline_filter(padded, order, output=numpy.float64, mode=mode) else: npad = 0 filtered = input mode = _ni_support._extend_mode_to_code(mode) matrix = numpy.asarray(matrix, dtype=numpy.float64) if matrix.ndim not in [1, 2] or matrix.shape[0] < 1: raise RuntimeError('no proper affine matrix provided') if (matrix.ndim == 2 and matrix.shape[1] == input.ndim + 1 and (matrix.shape[0] in [input.ndim, input.ndim + 1])): if matrix.shape[0] == input.ndim + 1: exptd = [0] * input.ndim + [1] if not numpy.all(matrix[input.ndim] == exptd): msg = ('Expected homogeneous transformation matrix with ' 'shape {} for image shape {}, but bottom row was ' 'not equal to {}'.format(matrix.shape, input.shape, exptd)) raise ValueError(msg) # assume input is homogeneous coordinate transformation matrix offset = matrix[:input.ndim, input.ndim] matrix = matrix[:input.ndim, :input.ndim] if matrix.shape[0] != input.ndim: raise RuntimeError('affine matrix has wrong number of rows') if matrix.ndim == 2 and matrix.shape[1] != output.ndim: raise RuntimeError('affine matrix has wrong number of columns') if not matrix.flags.contiguous: matrix = matrix.copy() offset = _ni_support._normalize_sequence(offset, input.ndim) offset = numpy.asarray(offset, dtype=numpy.float64) if offset.ndim != 1 or offset.shape[0] < 1: raise RuntimeError('no proper offset provided') if not offset.flags.contiguous: offset = offset.copy() if matrix.ndim == 1: warnings.warn( "The behavior of affine_transform with a 1-D " "array supplied for the matrix parameter has changed in " "SciPy 0.18.0.", stacklevel=2 ) _nd_image.zoom_shift(filtered, matrix, offset/matrix, output, order, mode, cval, npad, False) else: _nd_image.geometric_transform(filtered, None, None, matrix, offset, output, order, mode, cval, npad, None, None) return output @docfiller def shift(input, shift, output=None, order=3, mode='constant', cval=0.0, prefilter=True): """ Shift an array. The array is shifted using spline interpolation of the requested order. Points outside the boundaries of the input are filled according to the given mode. Parameters ---------- %(input)s shift : float or sequence The shift along the axes. If a float, `shift` is the same for each axis. If a sequence, `shift` should contain one value for each axis. %(output)s order : int, optional The order of the spline interpolation, default is 3. The order has to be in the range 0-5. %(mode_interp_constant)s %(cval)s %(prefilter)s Returns ------- shift : ndarray The shifted input. See Also -------- affine_transform : Affine transformations Notes ----- For complex-valued `input`, this function shifts the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- Import the necessary modules and an exemplary image. >>> from scipy.ndimage import shift >>> import matplotlib.pyplot as plt >>> from scipy import datasets >>> image = datasets.ascent() Shift the image vertically by 20 pixels. >>> image_shifted_vertically = shift(image, (20, 0)) Shift the image vertically by -200 pixels and horizontally by 100 pixels. >>> image_shifted_both_directions = shift(image, (-200, 100)) Plot the original and the shifted images. >>> fig, axes = plt.subplots(3, 1, figsize=(4, 12)) >>> plt.gray() # show the filtered result in grayscale >>> top, middle, bottom = axes >>> for ax in axes: ... ax.set_axis_off() # remove coordinate system >>> top.imshow(image) >>> top.set_title("Original image") >>> middle.imshow(image_shifted_vertically) >>> middle.set_title("Vertically shifted image") >>> bottom.imshow(image_shifted_both_directions) >>> bottom.set_title("Image shifted in both directions") >>> fig.tight_layout() """ if order < 0 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) if input.ndim < 1: raise RuntimeError('input and output rank must be > 0') complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, complex_output=complex_output) if complex_output: # import under different name to avoid confusion with shift parameter from scipy.ndimage._interpolation import shift as _shift kwargs = dict(order=order, mode=mode, prefilter=prefilter) _shift(input.real, shift, output=output.real, cval=numpy.real(cval), **kwargs) _shift(input.imag, shift, output=output.imag, cval=numpy.imag(cval), **kwargs) return output if prefilter and order > 1: padded, npad = _prepad_for_spline_filter(input, mode, cval) filtered = spline_filter(padded, order, output=numpy.float64, mode=mode) else: npad = 0 filtered = input mode = _ni_support._extend_mode_to_code(mode) shift = _ni_support._normalize_sequence(shift, input.ndim) shift = [-ii for ii in shift] shift = numpy.asarray(shift, dtype=numpy.float64) if not shift.flags.contiguous: shift = shift.copy() _nd_image.zoom_shift(filtered, None, shift, output, order, mode, cval, npad, False) return output @docfiller def zoom(input, zoom, output=None, order=3, mode='constant', cval=0.0, prefilter=True, *, grid_mode=False): """ Zoom an array. The array is zoomed using spline interpolation of the requested order. Parameters ---------- %(input)s zoom : float or sequence The zoom factor along the axes. If a float, `zoom` is the same for each axis. If a sequence, `zoom` should contain one value for each axis. %(output)s order : int, optional The order of the spline interpolation, default is 3. The order has to be in the range 0-5. %(mode_interp_constant)s %(cval)s %(prefilter)s grid_mode : bool, optional If False, the distance from the pixel centers is zoomed. Otherwise, the distance including the full pixel extent is used. For example, a 1d signal of length 5 is considered to have length 4 when `grid_mode` is False, but length 5 when `grid_mode` is True. See the following visual illustration: .. code-block:: text | pixel 1 | pixel 2 | pixel 3 | pixel 4 | pixel 5 | |<-------------------------------------->| vs. |<----------------------------------------------->| The starting point of the arrow in the diagram above corresponds to coordinate location 0 in each mode. Returns ------- zoom : ndarray The zoomed input. Notes ----- For complex-valued `input`, this function zooms the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- >>> from scipy import ndimage, datasets >>> import matplotlib.pyplot as plt >>> fig = plt.figure() >>> ax1 = fig.add_subplot(121) # left side >>> ax2 = fig.add_subplot(122) # right side >>> ascent = datasets.ascent() >>> result = ndimage.zoom(ascent, 3.0) >>> ax1.imshow(ascent, vmin=0, vmax=255) >>> ax2.imshow(result, vmin=0, vmax=255) >>> plt.show() >>> print(ascent.shape) (512, 512) >>> print(result.shape) (1536, 1536) """ if order < 0 or order > 5: raise RuntimeError('spline order not supported') input = numpy.asarray(input) if input.ndim < 1: raise RuntimeError('input and output rank must be > 0') zoom = _ni_support._normalize_sequence(zoom, input.ndim) output_shape = tuple( [int(round(ii * jj)) for ii, jj in zip(input.shape, zoom)]) complex_output = numpy.iscomplexobj(input) output = _ni_support._get_output(output, input, shape=output_shape, complex_output=complex_output) if complex_output: # import under different name to avoid confusion with zoom parameter from scipy.ndimage._interpolation import zoom as _zoom kwargs = dict(order=order, mode=mode, prefilter=prefilter) _zoom(input.real, zoom, output=output.real, cval=numpy.real(cval), **kwargs) _zoom(input.imag, zoom, output=output.imag, cval=numpy.imag(cval), **kwargs) return output if prefilter and order > 1: padded, npad = _prepad_for_spline_filter(input, mode, cval) filtered = spline_filter(padded, order, output=numpy.float64, mode=mode) else: npad = 0 filtered = input if grid_mode: # warn about modes that may have surprising behavior suggest_mode = None if mode == 'constant': suggest_mode = 'grid-constant' elif mode == 'wrap': suggest_mode = 'grid-wrap' if suggest_mode is not None: warnings.warn( ("It is recommended to use mode = {} instead of {} when " "grid_mode is True.").format(suggest_mode, mode), stacklevel=2 ) mode = _ni_support._extend_mode_to_code(mode) zoom_div = numpy.array(output_shape) zoom_nominator = numpy.array(input.shape) if not grid_mode: zoom_div -= 1 zoom_nominator -= 1 # Zooming to infinite values is unpredictable, so just choose # zoom factor 1 instead zoom = numpy.divide(zoom_nominator, zoom_div, out=numpy.ones_like(input.shape, dtype=numpy.float64), where=zoom_div != 0) zoom = numpy.ascontiguousarray(zoom) _nd_image.zoom_shift(filtered, zoom, None, output, order, mode, cval, npad, grid_mode) return output @docfiller def rotate(input, angle, axes=(1, 0), reshape=True, output=None, order=3, mode='constant', cval=0.0, prefilter=True): """ Rotate an array. The array is rotated in the plane defined by the two axes given by the `axes` parameter using spline interpolation of the requested order. Parameters ---------- %(input)s angle : float The rotation angle in degrees. axes : tuple of 2 ints, optional The two axes that define the plane of rotation. Default is the first two axes. reshape : bool, optional If `reshape` is true, the output shape is adapted so that the input array is contained completely in the output. Default is True. %(output)s order : int, optional The order of the spline interpolation, default is 3. The order has to be in the range 0-5. %(mode_interp_constant)s %(cval)s %(prefilter)s Returns ------- rotate : ndarray The rotated input. Notes ----- For complex-valued `input`, this function rotates the real and imaginary components independently. .. versionadded:: 1.6.0 Complex-valued support added. Examples -------- >>> from scipy import ndimage, datasets >>> import matplotlib.pyplot as plt >>> fig = plt.figure(figsize=(10, 3)) >>> ax1, ax2, ax3 = fig.subplots(1, 3) >>> img = datasets.ascent() >>> img_45 = ndimage.rotate(img, 45, reshape=False) >>> full_img_45 = ndimage.rotate(img, 45, reshape=True) >>> ax1.imshow(img, cmap='gray') >>> ax1.set_axis_off() >>> ax2.imshow(img_45, cmap='gray') >>> ax2.set_axis_off() >>> ax3.imshow(full_img_45, cmap='gray') >>> ax3.set_axis_off() >>> fig.set_layout_engine('tight') >>> plt.show() >>> print(img.shape) (512, 512) >>> print(img_45.shape) (512, 512) >>> print(full_img_45.shape) (724, 724) """ input_arr = numpy.asarray(input) ndim = input_arr.ndim if ndim < 2: raise ValueError('input array should be at least 2D') axes = list(axes) if len(axes) != 2: raise ValueError('axes should contain exactly two values') if not all([float(ax).is_integer() for ax in axes]): raise ValueError('axes should contain only integer values') if axes[0] < 0: axes[0] += ndim if axes[1] < 0: axes[1] += ndim if axes[0] < 0 or axes[1] < 0 or axes[0] >= ndim or axes[1] >= ndim: raise ValueError('invalid rotation plane specified') axes.sort() c, s = special.cosdg(angle), special.sindg(angle) rot_matrix = numpy.array([[c, s], [-s, c]]) img_shape = numpy.asarray(input_arr.shape) in_plane_shape = img_shape[axes] if reshape: # Compute transformed input bounds iy, ix = in_plane_shape out_bounds = rot_matrix @ [[0, 0, iy, iy], [0, ix, 0, ix]] # Compute the shape of the transformed input plane out_plane_shape = (numpy.ptp(out_bounds, axis=1) + 0.5).astype(int) else: out_plane_shape = img_shape[axes] out_center = rot_matrix @ ((out_plane_shape - 1) / 2) in_center = (in_plane_shape - 1) / 2 offset = in_center - out_center output_shape = img_shape output_shape[axes] = out_plane_shape output_shape = tuple(output_shape) complex_output = numpy.iscomplexobj(input_arr) output = _ni_support._get_output(output, input_arr, shape=output_shape, complex_output=complex_output) if ndim <= 2: affine_transform(input_arr, rot_matrix, offset, output_shape, output, order, mode, cval, prefilter) else: # If ndim > 2, the rotation is applied over all the planes # parallel to axes planes_coord = itertools.product( *[[slice(None)] if ax in axes else range(img_shape[ax]) for ax in range(ndim)]) out_plane_shape = tuple(out_plane_shape) for coordinates in planes_coord: ia = input_arr[coordinates] oa = output[coordinates] affine_transform(ia, rot_matrix, offset, out_plane_shape, oa, order, mode, cval, prefilter) return output