1"""
2Discrete Fourier Transforms
3
4Routines in this module:
5
6fft(a, n=None, axis=-1, norm="backward")
7ifft(a, n=None, axis=-1, norm="backward")
8rfft(a, n=None, axis=-1, norm="backward")
9irfft(a, n=None, axis=-1, norm="backward")
10hfft(a, n=None, axis=-1, norm="backward")
11ihfft(a, n=None, axis=-1, norm="backward")
12fftn(a, s=None, axes=None, norm="backward")
13ifftn(a, s=None, axes=None, norm="backward")
14rfftn(a, s=None, axes=None, norm="backward")
15irfftn(a, s=None, axes=None, norm="backward")
16fft2(a, s=None, axes=(-2,-1), norm="backward")
17ifft2(a, s=None, axes=(-2, -1), norm="backward")
18rfft2(a, s=None, axes=(-2,-1), norm="backward")
19irfft2(a, s=None, axes=(-2, -1), norm="backward")
20
21i = inverse transform
22r = transform of purely real data
23h = Hermite transform
24n = n-dimensional transform
252 = 2-dimensional transform
26(Note: 2D routines are just nD routines with different default
27behavior.)
28
29"""
30__all__ = ['fft', 'ifft', 'rfft', 'irfft', 'hfft', 'ihfft', 'rfftn',
31           'irfftn', 'rfft2', 'irfft2', 'fft2', 'ifft2', 'fftn', 'ifftn']
32
33import functools
34
35from numpy.core import asarray, zeros, swapaxes, conjugate, take, sqrt
36from . import _pocketfft_internal as pfi
37from numpy.core.multiarray import normalize_axis_index
38from numpy.core import overrides
39
40
41array_function_dispatch = functools.partial(
42    overrides.array_function_dispatch, module='numpy.fft')
43
44
45# `inv_norm` is a float by which the result of the transform needs to be
46# divided. This replaces the original, more intuitive 'fct` parameter to avoid
47# divisions by zero (or alternatively additional checks) in the case of
48# zero-length axes during its computation.
49def _raw_fft(a, n, axis, is_real, is_forward, inv_norm):
50    axis = normalize_axis_index(axis, a.ndim)
51    if n is None:
52        n = a.shape[axis]
53
54    fct = 1/inv_norm
55
56    if a.shape[axis] != n:
57        s = list(a.shape)
58        index = [slice(None)]*len(s)
59        if s[axis] > n:
60            index[axis] = slice(0, n)
61            a = a[tuple(index)]
62        else:
63            index[axis] = slice(0, s[axis])
64            s[axis] = n
65            z = zeros(s, a.dtype.char)
66            z[tuple(index)] = a
67            a = z
68
69    if axis == a.ndim-1:
70        r = pfi.execute(a, is_real, is_forward, fct)
71    else:
72        a = swapaxes(a, axis, -1)
73        r = pfi.execute(a, is_real, is_forward, fct)
74        r = swapaxes(r, axis, -1)
75    return r
76
77
78def _get_forward_norm(n, norm):
79    if n < 1:
80        raise ValueError(f"Invalid number of FFT data points ({n}) specified.")
81
82    if norm is None or norm == "backward":
83        return 1
84    elif norm == "ortho":
85        return sqrt(n)
86    elif norm == "forward":
87        return n
88    raise ValueError(f'Invalid norm value {norm}; should be "backward",'
89                     '"ortho" or "forward".')
90
91
92def _get_backward_norm(n, norm):
93    if n < 1:
94        raise ValueError(f"Invalid number of FFT data points ({n}) specified.")
95
96    if norm is None or norm == "backward":
97        return n
98    elif norm == "ortho":
99        return sqrt(n)
100    elif norm == "forward":
101        return 1
102    raise ValueError(f'Invalid norm value {norm}; should be "backward", '
103                     '"ortho" or "forward".')
104
105
106_SWAP_DIRECTION_MAP = {"backward": "forward", None: "forward",
107                       "ortho": "ortho", "forward": "backward"}
108
109
110def _swap_direction(norm):
111    try:
112        return _SWAP_DIRECTION_MAP[norm]
113    except KeyError:
114        raise ValueError(f'Invalid norm value {norm}; should be "backward", '
115                         '"ortho" or "forward".')
116
117
118def _fft_dispatcher(a, n=None, axis=None, norm=None):
119    return (a,)
120
121
122@array_function_dispatch(_fft_dispatcher)
123def fft(a, n=None, axis=-1, norm=None):
124    """
125    Compute the one-dimensional discrete Fourier Transform.
126
127    This function computes the one-dimensional *n*-point discrete Fourier
128    Transform (DFT) with the efficient Fast Fourier Transform (FFT)
129    algorithm [CT].
130
131    Parameters
132    ----------
133    a : array_like
134        Input array, can be complex.
135    n : int, optional
136        Length of the transformed axis of the output.
137        If `n` is smaller than the length of the input, the input is cropped.
138        If it is larger, the input is padded with zeros.  If `n` is not given,
139        the length of the input along the axis specified by `axis` is used.
140    axis : int, optional
141        Axis over which to compute the FFT.  If not given, the last axis is
142        used.
143    norm : {"backward", "ortho", "forward"}, optional
144        .. versionadded:: 1.10.0
145
146        Normalization mode (see `numpy.fft`). Default is "backward".
147        Indicates which direction of the forward/backward pair of transforms
148        is scaled and with what normalization factor.
149
150        .. versionadded:: 1.20.0
151
152            The "backward", "forward" values were added.
153
154    Returns
155    -------
156    out : complex ndarray
157        The truncated or zero-padded input, transformed along the axis
158        indicated by `axis`, or the last one if `axis` is not specified.
159
160    Raises
161    ------
162    IndexError
163        if `axes` is larger than the last axis of `a`.
164
165    See Also
166    --------
167    numpy.fft : for definition of the DFT and conventions used.
168    ifft : The inverse of `fft`.
169    fft2 : The two-dimensional FFT.
170    fftn : The *n*-dimensional FFT.
171    rfftn : The *n*-dimensional FFT of real input.
172    fftfreq : Frequency bins for given FFT parameters.
173
174    Notes
175    -----
176    FFT (Fast Fourier Transform) refers to a way the discrete Fourier
177    Transform (DFT) can be calculated efficiently, by using symmetries in the
178    calculated terms.  The symmetry is highest when `n` is a power of 2, and
179    the transform is therefore most efficient for these sizes.
180
181    The DFT is defined, with the conventions used in this implementation, in
182    the documentation for the `numpy.fft` module.
183
184    References
185    ----------
186    .. [CT] Cooley, James W., and John W. Tukey, 1965, "An algorithm for the
187            machine calculation of complex Fourier series," *Math. Comput.*
188            19: 297-301.
189
190    Examples
191    --------
192    >>> np.fft.fft(np.exp(2j * np.pi * np.arange(8) / 8))
193    array([-2.33486982e-16+1.14423775e-17j,  8.00000000e+00-1.25557246e-15j,
194            2.33486982e-16+2.33486982e-16j,  0.00000000e+00+1.22464680e-16j,
195           -1.14423775e-17+2.33486982e-16j,  0.00000000e+00+5.20784380e-16j,
196            1.14423775e-17+1.14423775e-17j,  0.00000000e+00+1.22464680e-16j])
197
198    In this example, real input has an FFT which is Hermitian, i.e., symmetric
199    in the real part and anti-symmetric in the imaginary part, as described in
200    the `numpy.fft` documentation:
201
202    >>> import matplotlib.pyplot as plt
203    >>> t = np.arange(256)
204    >>> sp = np.fft.fft(np.sin(t))
205    >>> freq = np.fft.fftfreq(t.shape[-1])
206    >>> plt.plot(freq, sp.real, freq, sp.imag)
207    [<matplotlib.lines.Line2D object at 0x...>, <matplotlib.lines.Line2D object at 0x...>]
208    >>> plt.show()
209
210    """
211    a = asarray(a)
212    if n is None:
213        n = a.shape[axis]
214    inv_norm = _get_forward_norm(n, norm)
215    output = _raw_fft(a, n, axis, False, True, inv_norm)
216    return output
217
218
219@array_function_dispatch(_fft_dispatcher)
220def ifft(a, n=None, axis=-1, norm=None):
221    """
222    Compute the one-dimensional inverse discrete Fourier Transform.
223
224    This function computes the inverse of the one-dimensional *n*-point
225    discrete Fourier transform computed by `fft`.  In other words,
226    ``ifft(fft(a)) == a`` to within numerical accuracy.
227    For a general description of the algorithm and definitions,
228    see `numpy.fft`.
229
230    The input should be ordered in the same way as is returned by `fft`,
231    i.e.,
232
233    * ``a[0]`` should contain the zero frequency term,
234    * ``a[1:n//2]`` should contain the positive-frequency terms,
235    * ``a[n//2 + 1:]`` should contain the negative-frequency terms, in
236      increasing order starting from the most negative frequency.
237
238    For an even number of input points, ``A[n//2]`` represents the sum of
239    the values at the positive and negative Nyquist frequencies, as the two
240    are aliased together. See `numpy.fft` for details.
241
242    Parameters
243    ----------
244    a : array_like
245        Input array, can be complex.
246    n : int, optional
247        Length of the transformed axis of the output.
248        If `n` is smaller than the length of the input, the input is cropped.
249        If it is larger, the input is padded with zeros.  If `n` is not given,
250        the length of the input along the axis specified by `axis` is used.
251        See notes about padding issues.
252    axis : int, optional
253        Axis over which to compute the inverse DFT.  If not given, the last
254        axis is used.
255    norm : {"backward", "ortho", "forward"}, optional
256        .. versionadded:: 1.10.0
257
258        Normalization mode (see `numpy.fft`). Default is "backward".
259        Indicates which direction of the forward/backward pair of transforms
260        is scaled and with what normalization factor.
261
262        .. versionadded:: 1.20.0
263
264            The "backward", "forward" values were added.
265
266    Returns
267    -------
268    out : complex ndarray
269        The truncated or zero-padded input, transformed along the axis
270        indicated by `axis`, or the last one if `axis` is not specified.
271
272    Raises
273    ------
274    IndexError
275        If `axes` is larger than the last axis of `a`.
276
277    See Also
278    --------
279    numpy.fft : An introduction, with definitions and general explanations.
280    fft : The one-dimensional (forward) FFT, of which `ifft` is the inverse
281    ifft2 : The two-dimensional inverse FFT.
282    ifftn : The n-dimensional inverse FFT.
283
284    Notes
285    -----
286    If the input parameter `n` is larger than the size of the input, the input
287    is padded by appending zeros at the end.  Even though this is the common
288    approach, it might lead to surprising results.  If a different padding is
289    desired, it must be performed before calling `ifft`.
290
291    Examples
292    --------
293    >>> np.fft.ifft([0, 4, 0, 0])
294    array([ 1.+0.j,  0.+1.j, -1.+0.j,  0.-1.j]) # may vary
295
296    Create and plot a band-limited signal with random phases:
297
298    >>> import matplotlib.pyplot as plt
299    >>> t = np.arange(400)
300    >>> n = np.zeros((400,), dtype=complex)
301    >>> n[40:60] = np.exp(1j*np.random.uniform(0, 2*np.pi, (20,)))
302    >>> s = np.fft.ifft(n)
303    >>> plt.plot(t, s.real, 'b-', t, s.imag, 'r--')
304    [<matplotlib.lines.Line2D object at ...>, <matplotlib.lines.Line2D object at ...>]
305    >>> plt.legend(('real', 'imaginary'))
306    <matplotlib.legend.Legend object at ...>
307    >>> plt.show()
308
309    """
310    a = asarray(a)
311    if n is None:
312        n = a.shape[axis]
313    inv_norm = _get_backward_norm(n, norm)
314    output = _raw_fft(a, n, axis, False, False, inv_norm)
315    return output
316
317
318@array_function_dispatch(_fft_dispatcher)
319def rfft(a, n=None, axis=-1, norm=None):
320    """
321    Compute the one-dimensional discrete Fourier Transform for real input.
322
323    This function computes the one-dimensional *n*-point discrete Fourier
324    Transform (DFT) of a real-valued array by means of an efficient algorithm
325    called the Fast Fourier Transform (FFT).
326
327    Parameters
328    ----------
329    a : array_like
330        Input array
331    n : int, optional
332        Number of points along transformation axis in the input to use.
333        If `n` is smaller than the length of the input, the input is cropped.
334        If it is larger, the input is padded with zeros. If `n` is not given,
335        the length of the input along the axis specified by `axis` is used.
336    axis : int, optional
337        Axis over which to compute the FFT. If not given, the last axis is
338        used.
339    norm : {"backward", "ortho", "forward"}, optional
340        .. versionadded:: 1.10.0
341
342        Normalization mode (see `numpy.fft`). Default is "backward".
343        Indicates which direction of the forward/backward pair of transforms
344        is scaled and with what normalization factor.
345
346        .. versionadded:: 1.20.0
347
348            The "backward", "forward" values were added.
349
350    Returns
351    -------
352    out : complex ndarray
353        The truncated or zero-padded input, transformed along the axis
354        indicated by `axis`, or the last one if `axis` is not specified.
355        If `n` is even, the length of the transformed axis is ``(n/2)+1``.
356        If `n` is odd, the length is ``(n+1)/2``.
357
358    Raises
359    ------
360    IndexError
361        If `axis` is larger than the last axis of `a`.
362
363    See Also
364    --------
365    numpy.fft : For definition of the DFT and conventions used.
366    irfft : The inverse of `rfft`.
367    fft : The one-dimensional FFT of general (complex) input.
368    fftn : The *n*-dimensional FFT.
369    rfftn : The *n*-dimensional FFT of real input.
370
371    Notes
372    -----
373    When the DFT is computed for purely real input, the output is
374    Hermitian-symmetric, i.e. the negative frequency terms are just the complex
375    conjugates of the corresponding positive-frequency terms, and the
376    negative-frequency terms are therefore redundant.  This function does not
377    compute the negative frequency terms, and the length of the transformed
378    axis of the output is therefore ``n//2 + 1``.
379
380    When ``A = rfft(a)`` and fs is the sampling frequency, ``A[0]`` contains
381    the zero-frequency term 0*fs, which is real due to Hermitian symmetry.
382
383    If `n` is even, ``A[-1]`` contains the term representing both positive
384    and negative Nyquist frequency (+fs/2 and -fs/2), and must also be purely
385    real. If `n` is odd, there is no term at fs/2; ``A[-1]`` contains
386    the largest positive frequency (fs/2*(n-1)/n), and is complex in the
387    general case.
388
389    If the input `a` contains an imaginary part, it is silently discarded.
390
391    Examples
392    --------
393    >>> np.fft.fft([0, 1, 0, 0])
394    array([ 1.+0.j,  0.-1.j, -1.+0.j,  0.+1.j]) # may vary
395    >>> np.fft.rfft([0, 1, 0, 0])
396    array([ 1.+0.j,  0.-1.j, -1.+0.j]) # may vary
397
398    Notice how the final element of the `fft` output is the complex conjugate
399    of the second element, for real input. For `rfft`, this symmetry is
400    exploited to compute only the non-negative frequency terms.
401
402    """
403    a = asarray(a)
404    if n is None:
405        n = a.shape[axis]
406    inv_norm = _get_forward_norm(n, norm)
407    output = _raw_fft(a, n, axis, True, True, inv_norm)
408    return output
409
410
411@array_function_dispatch(_fft_dispatcher)
412def irfft(a, n=None, axis=-1, norm=None):
413    """
414    Computes the inverse of `rfft`.
415
416    This function computes the inverse of the one-dimensional *n*-point
417    discrete Fourier Transform of real input computed by `rfft`.
418    In other words, ``irfft(rfft(a), len(a)) == a`` to within numerical
419    accuracy. (See Notes below for why ``len(a)`` is necessary here.)
420
421    The input is expected to be in the form returned by `rfft`, i.e. the
422    real zero-frequency term followed by the complex positive frequency terms
423    in order of increasing frequency.  Since the discrete Fourier Transform of
424    real input is Hermitian-symmetric, the negative frequency terms are taken
425    to be the complex conjugates of the corresponding positive frequency terms.
426
427    Parameters
428    ----------
429    a : array_like
430        The input array.
431    n : int, optional
432        Length of the transformed axis of the output.
433        For `n` output points, ``n//2+1`` input points are necessary.  If the
434        input is longer than this, it is cropped.  If it is shorter than this,
435        it is padded with zeros.  If `n` is not given, it is taken to be
436        ``2*(m-1)`` where ``m`` is the length of the input along the axis
437        specified by `axis`.
438    axis : int, optional
439        Axis over which to compute the inverse FFT. If not given, the last
440        axis is used.
441    norm : {"backward", "ortho", "forward"}, optional
442        .. versionadded:: 1.10.0
443
444        Normalization mode (see `numpy.fft`). Default is "backward".
445        Indicates which direction of the forward/backward pair of transforms
446        is scaled and with what normalization factor.
447
448        .. versionadded:: 1.20.0
449
450            The "backward", "forward" values were added.
451
452    Returns
453    -------
454    out : ndarray
455        The truncated or zero-padded input, transformed along the axis
456        indicated by `axis`, or the last one if `axis` is not specified.
457        The length of the transformed axis is `n`, or, if `n` is not given,
458        ``2*(m-1)`` where ``m`` is the length of the transformed axis of the
459        input. To get an odd number of output points, `n` must be specified.
460
461    Raises
462    ------
463    IndexError
464        If `axis` is larger than the last axis of `a`.
465
466    See Also
467    --------
468    numpy.fft : For definition of the DFT and conventions used.
469    rfft : The one-dimensional FFT of real input, of which `irfft` is inverse.
470    fft : The one-dimensional FFT.
471    irfft2 : The inverse of the two-dimensional FFT of real input.
472    irfftn : The inverse of the *n*-dimensional FFT of real input.
473
474    Notes
475    -----
476    Returns the real valued `n`-point inverse discrete Fourier transform
477    of `a`, where `a` contains the non-negative frequency terms of a
478    Hermitian-symmetric sequence. `n` is the length of the result, not the
479    input.
480
481    If you specify an `n` such that `a` must be zero-padded or truncated, the
482    extra/removed values will be added/removed at high frequencies. One can
483    thus resample a series to `m` points via Fourier interpolation by:
484    ``a_resamp = irfft(rfft(a), m)``.
485
486    The correct interpretation of the hermitian input depends on the length of
487    the original data, as given by `n`. This is because each input shape could
488    correspond to either an odd or even length signal. By default, `irfft`
489    assumes an even output length which puts the last entry at the Nyquist
490    frequency; aliasing with its symmetric counterpart. By Hermitian symmetry,
491    the value is thus treated as purely real. To avoid losing information, the
492    correct length of the real input **must** be given.
493
494    Examples
495    --------
496    >>> np.fft.ifft([1, -1j, -1, 1j])
497    array([0.+0.j,  1.+0.j,  0.+0.j,  0.+0.j]) # may vary
498    >>> np.fft.irfft([1, -1j, -1])
499    array([0.,  1.,  0.,  0.])
500
501    Notice how the last term in the input to the ordinary `ifft` is the
502    complex conjugate of the second term, and the output has zero imaginary
503    part everywhere.  When calling `irfft`, the negative frequencies are not
504    specified, and the output array is purely real.
505
506    """
507    a = asarray(a)
508    if n is None:
509        n = (a.shape[axis] - 1) * 2
510    inv_norm = _get_backward_norm(n, norm)
511    output = _raw_fft(a, n, axis, True, False, inv_norm)
512    return output
513
514
515@array_function_dispatch(_fft_dispatcher)
516def hfft(a, n=None, axis=-1, norm=None):
517    """
518    Compute the FFT of a signal that has Hermitian symmetry, i.e., a real
519    spectrum.
520
521    Parameters
522    ----------
523    a : array_like
524        The input array.
525    n : int, optional
526        Length of the transformed axis of the output. For `n` output
527        points, ``n//2 + 1`` input points are necessary.  If the input is
528        longer than this, it is cropped.  If it is shorter than this, it is
529        padded with zeros.  If `n` is not given, it is taken to be ``2*(m-1)``
530        where ``m`` is the length of the input along the axis specified by
531        `axis`.
532    axis : int, optional
533        Axis over which to compute the FFT. If not given, the last
534        axis is used.
535    norm : {"backward", "ortho", "forward"}, optional
536        .. versionadded:: 1.10.0
537
538        Normalization mode (see `numpy.fft`). Default is "backward".
539        Indicates which direction of the forward/backward pair of transforms
540        is scaled and with what normalization factor.
541
542        .. versionadded:: 1.20.0
543
544            The "backward", "forward" values were added.
545
546    Returns
547    -------
548    out : ndarray
549        The truncated or zero-padded input, transformed along the axis
550        indicated by `axis`, or the last one if `axis` is not specified.
551        The length of the transformed axis is `n`, or, if `n` is not given,
552        ``2*m - 2`` where ``m`` is the length of the transformed axis of
553        the input. To get an odd number of output points, `n` must be
554        specified, for instance as ``2*m - 1`` in the typical case,
555
556    Raises
557    ------
558    IndexError
559        If `axis` is larger than the last axis of `a`.
560
561    See also
562    --------
563    rfft : Compute the one-dimensional FFT for real input.
564    ihfft : The inverse of `hfft`.
565
566    Notes
567    -----
568    `hfft`/`ihfft` are a pair analogous to `rfft`/`irfft`, but for the
569    opposite case: here the signal has Hermitian symmetry in the time
570    domain and is real in the frequency domain. So here it's `hfft` for
571    which you must supply the length of the result if it is to be odd.
572
573    * even: ``ihfft(hfft(a, 2*len(a) - 2)) == a``, within roundoff error,
574    * odd: ``ihfft(hfft(a, 2*len(a) - 1)) == a``, within roundoff error.
575
576    The correct interpretation of the hermitian input depends on the length of
577    the original data, as given by `n`. This is because each input shape could
578    correspond to either an odd or even length signal. By default, `hfft`
579    assumes an even output length which puts the last entry at the Nyquist
580    frequency; aliasing with its symmetric counterpart. By Hermitian symmetry,
581    the value is thus treated as purely real. To avoid losing information, the
582    shape of the full signal **must** be given.
583
584    Examples
585    --------
586    >>> signal = np.array([1, 2, 3, 4, 3, 2])
587    >>> np.fft.fft(signal)
588    array([15.+0.j,  -4.+0.j,   0.+0.j,  -1.-0.j,   0.+0.j,  -4.+0.j]) # may vary
589    >>> np.fft.hfft(signal[:4]) # Input first half of signal
590    array([15.,  -4.,   0.,  -1.,   0.,  -4.])
591    >>> np.fft.hfft(signal, 6)  # Input entire signal and truncate
592    array([15.,  -4.,   0.,  -1.,   0.,  -4.])
593
594
595    >>> signal = np.array([[1, 1.j], [-1.j, 2]])
596    >>> np.conj(signal.T) - signal   # check Hermitian symmetry
597    array([[ 0.-0.j,  -0.+0.j], # may vary
598           [ 0.+0.j,  0.-0.j]])
599    >>> freq_spectrum = np.fft.hfft(signal)
600    >>> freq_spectrum
601    array([[ 1.,  1.],
602           [ 2., -2.]])
603
604    """
605    a = asarray(a)
606    if n is None:
607        n = (a.shape[axis] - 1) * 2
608    new_norm = _swap_direction(norm)
609    output = irfft(conjugate(a), n, axis, norm=new_norm)
610    return output
611
612
613@array_function_dispatch(_fft_dispatcher)
614def ihfft(a, n=None, axis=-1, norm=None):
615    """
616    Compute the inverse FFT of a signal that has Hermitian symmetry.
617
618    Parameters
619    ----------
620    a : array_like
621        Input array.
622    n : int, optional
623        Length of the inverse FFT, the number of points along
624        transformation axis in the input to use.  If `n` is smaller than
625        the length of the input, the input is cropped.  If it is larger,
626        the input is padded with zeros. If `n` is not given, the length of
627        the input along the axis specified by `axis` is used.
628    axis : int, optional
629        Axis over which to compute the inverse FFT. If not given, the last
630        axis is used.
631    norm : {"backward", "ortho", "forward"}, optional
632        .. versionadded:: 1.10.0
633
634        Normalization mode (see `numpy.fft`). Default is "backward".
635        Indicates which direction of the forward/backward pair of transforms
636        is scaled and with what normalization factor.
637
638        .. versionadded:: 1.20.0
639
640            The "backward", "forward" values were added.
641
642    Returns
643    -------
644    out : complex ndarray
645        The truncated or zero-padded input, transformed along the axis
646        indicated by `axis`, or the last one if `axis` is not specified.
647        The length of the transformed axis is ``n//2 + 1``.
648
649    See also
650    --------
651    hfft, irfft
652
653    Notes
654    -----
655    `hfft`/`ihfft` are a pair analogous to `rfft`/`irfft`, but for the
656    opposite case: here the signal has Hermitian symmetry in the time
657    domain and is real in the frequency domain. So here it's `hfft` for
658    which you must supply the length of the result if it is to be odd:
659
660    * even: ``ihfft(hfft(a, 2*len(a) - 2)) == a``, within roundoff error,
661    * odd: ``ihfft(hfft(a, 2*len(a) - 1)) == a``, within roundoff error.
662
663    Examples
664    --------
665    >>> spectrum = np.array([ 15, -4, 0, -1, 0, -4])
666    >>> np.fft.ifft(spectrum)
667    array([1.+0.j,  2.+0.j,  3.+0.j,  4.+0.j,  3.+0.j,  2.+0.j]) # may vary
668    >>> np.fft.ihfft(spectrum)
669    array([ 1.-0.j,  2.-0.j,  3.-0.j,  4.-0.j]) # may vary
670
671    """
672    a = asarray(a)
673    if n is None:
674        n = a.shape[axis]
675    new_norm = _swap_direction(norm)
676    output = conjugate(rfft(a, n, axis, norm=new_norm))
677    return output
678
679
680def _cook_nd_args(a, s=None, axes=None, invreal=0):
681    if s is None:
682        shapeless = 1
683        if axes is None:
684            s = list(a.shape)
685        else:
686            s = take(a.shape, axes)
687    else:
688        shapeless = 0
689    s = list(s)
690    if axes is None:
691        axes = list(range(-len(s), 0))
692    if len(s) != len(axes):
693        raise ValueError("Shape and axes have different lengths.")
694    if invreal and shapeless:
695        s[-1] = (a.shape[axes[-1]] - 1) * 2
696    return s, axes
697
698
699def _raw_fftnd(a, s=None, axes=None, function=fft, norm=None):
700    a = asarray(a)
701    s, axes = _cook_nd_args(a, s, axes)
702    itl = list(range(len(axes)))
703    itl.reverse()
704    for ii in itl:
705        a = function(a, n=s[ii], axis=axes[ii], norm=norm)
706    return a
707
708
709def _fftn_dispatcher(a, s=None, axes=None, norm=None):
710    return (a,)
711
712
713@array_function_dispatch(_fftn_dispatcher)
714def fftn(a, s=None, axes=None, norm=None):
715    """
716    Compute the N-dimensional discrete Fourier Transform.
717
718    This function computes the *N*-dimensional discrete Fourier Transform over
719    any number of axes in an *M*-dimensional array by means of the Fast Fourier
720    Transform (FFT).
721
722    Parameters
723    ----------
724    a : array_like
725        Input array, can be complex.
726    s : sequence of ints, optional
727        Shape (length of each transformed axis) of the output
728        (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
729        This corresponds to ``n`` for ``fft(x, n)``.
730        Along any axis, if the given shape is smaller than that of the input,
731        the input is cropped.  If it is larger, the input is padded with zeros.
732        if `s` is not given, the shape of the input along the axes specified
733        by `axes` is used.
734    axes : sequence of ints, optional
735        Axes over which to compute the FFT.  If not given, the last ``len(s)``
736        axes are used, or all axes if `s` is also not specified.
737        Repeated indices in `axes` means that the transform over that axis is
738        performed multiple times.
739    norm : {"backward", "ortho", "forward"}, optional
740        .. versionadded:: 1.10.0
741
742        Normalization mode (see `numpy.fft`). Default is "backward".
743        Indicates which direction of the forward/backward pair of transforms
744        is scaled and with what normalization factor.
745
746        .. versionadded:: 1.20.0
747
748            The "backward", "forward" values were added.
749
750    Returns
751    -------
752    out : complex ndarray
753        The truncated or zero-padded input, transformed along the axes
754        indicated by `axes`, or by a combination of `s` and `a`,
755        as explained in the parameters section above.
756
757    Raises
758    ------
759    ValueError
760        If `s` and `axes` have different length.
761    IndexError
762        If an element of `axes` is larger than than the number of axes of `a`.
763
764    See Also
765    --------
766    numpy.fft : Overall view of discrete Fourier transforms, with definitions
767        and conventions used.
768    ifftn : The inverse of `fftn`, the inverse *n*-dimensional FFT.
769    fft : The one-dimensional FFT, with definitions and conventions used.
770    rfftn : The *n*-dimensional FFT of real input.
771    fft2 : The two-dimensional FFT.
772    fftshift : Shifts zero-frequency terms to centre of array
773
774    Notes
775    -----
776    The output, analogously to `fft`, contains the term for zero frequency in
777    the low-order corner of all axes, the positive frequency terms in the
778    first half of all axes, the term for the Nyquist frequency in the middle
779    of all axes and the negative frequency terms in the second half of all
780    axes, in order of decreasingly negative frequency.
781
782    See `numpy.fft` for details, definitions and conventions used.
783
784    Examples
785    --------
786    >>> a = np.mgrid[:3, :3, :3][0]
787    >>> np.fft.fftn(a, axes=(1, 2))
788    array([[[ 0.+0.j,   0.+0.j,   0.+0.j], # may vary
789            [ 0.+0.j,   0.+0.j,   0.+0.j],
790            [ 0.+0.j,   0.+0.j,   0.+0.j]],
791           [[ 9.+0.j,   0.+0.j,   0.+0.j],
792            [ 0.+0.j,   0.+0.j,   0.+0.j],
793            [ 0.+0.j,   0.+0.j,   0.+0.j]],
794           [[18.+0.j,   0.+0.j,   0.+0.j],
795            [ 0.+0.j,   0.+0.j,   0.+0.j],
796            [ 0.+0.j,   0.+0.j,   0.+0.j]]])
797    >>> np.fft.fftn(a, (2, 2), axes=(0, 1))
798    array([[[ 2.+0.j,  2.+0.j,  2.+0.j], # may vary
799            [ 0.+0.j,  0.+0.j,  0.+0.j]],
800           [[-2.+0.j, -2.+0.j, -2.+0.j],
801            [ 0.+0.j,  0.+0.j,  0.+0.j]]])
802
803    >>> import matplotlib.pyplot as plt
804    >>> [X, Y] = np.meshgrid(2 * np.pi * np.arange(200) / 12,
805    ...                      2 * np.pi * np.arange(200) / 34)
806    >>> S = np.sin(X) + np.cos(Y) + np.random.uniform(0, 1, X.shape)
807    >>> FS = np.fft.fftn(S)
808    >>> plt.imshow(np.log(np.abs(np.fft.fftshift(FS))**2))
809    <matplotlib.image.AxesImage object at 0x...>
810    >>> plt.show()
811
812    """
813    return _raw_fftnd(a, s, axes, fft, norm)
814
815
816@array_function_dispatch(_fftn_dispatcher)
817def ifftn(a, s=None, axes=None, norm=None):
818    """
819    Compute the N-dimensional inverse discrete Fourier Transform.
820
821    This function computes the inverse of the N-dimensional discrete
822    Fourier Transform over any number of axes in an M-dimensional array by
823    means of the Fast Fourier Transform (FFT).  In other words,
824    ``ifftn(fftn(a)) == a`` to within numerical accuracy.
825    For a description of the definitions and conventions used, see `numpy.fft`.
826
827    The input, analogously to `ifft`, should be ordered in the same way as is
828    returned by `fftn`, i.e. it should have the term for zero frequency
829    in all axes in the low-order corner, the positive frequency terms in the
830    first half of all axes, the term for the Nyquist frequency in the middle
831    of all axes and the negative frequency terms in the second half of all
832    axes, in order of decreasingly negative frequency.
833
834    Parameters
835    ----------
836    a : array_like
837        Input array, can be complex.
838    s : sequence of ints, optional
839        Shape (length of each transformed axis) of the output
840        (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
841        This corresponds to ``n`` for ``ifft(x, n)``.
842        Along any axis, if the given shape is smaller than that of the input,
843        the input is cropped.  If it is larger, the input is padded with zeros.
844        if `s` is not given, the shape of the input along the axes specified
845        by `axes` is used.  See notes for issue on `ifft` zero padding.
846    axes : sequence of ints, optional
847        Axes over which to compute the IFFT.  If not given, the last ``len(s)``
848        axes are used, or all axes if `s` is also not specified.
849        Repeated indices in `axes` means that the inverse transform over that
850        axis is performed multiple times.
851    norm : {"backward", "ortho", "forward"}, optional
852        .. versionadded:: 1.10.0
853
854        Normalization mode (see `numpy.fft`). Default is "backward".
855        Indicates which direction of the forward/backward pair of transforms
856        is scaled and with what normalization factor.
857
858        .. versionadded:: 1.20.0
859
860            The "backward", "forward" values were added.
861
862    Returns
863    -------
864    out : complex ndarray
865        The truncated or zero-padded input, transformed along the axes
866        indicated by `axes`, or by a combination of `s` or `a`,
867        as explained in the parameters section above.
868
869    Raises
870    ------
871    ValueError
872        If `s` and `axes` have different length.
873    IndexError
874        If an element of `axes` is larger than than the number of axes of `a`.
875
876    See Also
877    --------
878    numpy.fft : Overall view of discrete Fourier transforms, with definitions
879         and conventions used.
880    fftn : The forward *n*-dimensional FFT, of which `ifftn` is the inverse.
881    ifft : The one-dimensional inverse FFT.
882    ifft2 : The two-dimensional inverse FFT.
883    ifftshift : Undoes `fftshift`, shifts zero-frequency terms to beginning
884        of array.
885
886    Notes
887    -----
888    See `numpy.fft` for definitions and conventions used.
889
890    Zero-padding, analogously with `ifft`, is performed by appending zeros to
891    the input along the specified dimension.  Although this is the common
892    approach, it might lead to surprising results.  If another form of zero
893    padding is desired, it must be performed before `ifftn` is called.
894
895    Examples
896    --------
897    >>> a = np.eye(4)
898    >>> np.fft.ifftn(np.fft.fftn(a, axes=(0,)), axes=(1,))
899    array([[1.+0.j,  0.+0.j,  0.+0.j,  0.+0.j], # may vary
900           [0.+0.j,  1.+0.j,  0.+0.j,  0.+0.j],
901           [0.+0.j,  0.+0.j,  1.+0.j,  0.+0.j],
902           [0.+0.j,  0.+0.j,  0.+0.j,  1.+0.j]])
903
904
905    Create and plot an image with band-limited frequency content:
906
907    >>> import matplotlib.pyplot as plt
908    >>> n = np.zeros((200,200), dtype=complex)
909    >>> n[60:80, 20:40] = np.exp(1j*np.random.uniform(0, 2*np.pi, (20, 20)))
910    >>> im = np.fft.ifftn(n).real
911    >>> plt.imshow(im)
912    <matplotlib.image.AxesImage object at 0x...>
913    >>> plt.show()
914
915    """
916    return _raw_fftnd(a, s, axes, ifft, norm)
917
918
919@array_function_dispatch(_fftn_dispatcher)
920def fft2(a, s=None, axes=(-2, -1), norm=None):
921    """
922    Compute the 2-dimensional discrete Fourier Transform.
923
924    This function computes the *n*-dimensional discrete Fourier Transform
925    over any axes in an *M*-dimensional array by means of the
926    Fast Fourier Transform (FFT).  By default, the transform is computed over
927    the last two axes of the input array, i.e., a 2-dimensional FFT.
928
929    Parameters
930    ----------
931    a : array_like
932        Input array, can be complex
933    s : sequence of ints, optional
934        Shape (length of each transformed axis) of the output
935        (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
936        This corresponds to ``n`` for ``fft(x, n)``.
937        Along each axis, if the given shape is smaller than that of the input,
938        the input is cropped.  If it is larger, the input is padded with zeros.
939        if `s` is not given, the shape of the input along the axes specified
940        by `axes` is used.
941    axes : sequence of ints, optional
942        Axes over which to compute the FFT.  If not given, the last two
943        axes are used.  A repeated index in `axes` means the transform over
944        that axis is performed multiple times.  A one-element sequence means
945        that a one-dimensional FFT is performed.
946    norm : {"backward", "ortho", "forward"}, optional
947        .. versionadded:: 1.10.0
948
949        Normalization mode (see `numpy.fft`). Default is "backward".
950        Indicates which direction of the forward/backward pair of transforms
951        is scaled and with what normalization factor.
952
953        .. versionadded:: 1.20.0
954
955            The "backward", "forward" values were added.
956
957    Returns
958    -------
959    out : complex ndarray
960        The truncated or zero-padded input, transformed along the axes
961        indicated by `axes`, or the last two axes if `axes` is not given.
962
963    Raises
964    ------
965    ValueError
966        If `s` and `axes` have different length, or `axes` not given and
967        ``len(s) != 2``.
968    IndexError
969        If an element of `axes` is larger than than the number of axes of `a`.
970
971    See Also
972    --------
973    numpy.fft : Overall view of discrete Fourier transforms, with definitions
974         and conventions used.
975    ifft2 : The inverse two-dimensional FFT.
976    fft : The one-dimensional FFT.
977    fftn : The *n*-dimensional FFT.
978    fftshift : Shifts zero-frequency terms to the center of the array.
979        For two-dimensional input, swaps first and third quadrants, and second
980        and fourth quadrants.
981
982    Notes
983    -----
984    `fft2` is just `fftn` with a different default for `axes`.
985
986    The output, analogously to `fft`, contains the term for zero frequency in
987    the low-order corner of the transformed axes, the positive frequency terms
988    in the first half of these axes, the term for the Nyquist frequency in the
989    middle of the axes and the negative frequency terms in the second half of
990    the axes, in order of decreasingly negative frequency.
991
992    See `fftn` for details and a plotting example, and `numpy.fft` for
993    definitions and conventions used.
994
995
996    Examples
997    --------
998    >>> a = np.mgrid[:5, :5][0]
999    >>> np.fft.fft2(a)
1000    array([[ 50.  +0.j        ,   0.  +0.j        ,   0.  +0.j        , # may vary
1001              0.  +0.j        ,   0.  +0.j        ],
1002           [-12.5+17.20477401j,   0.  +0.j        ,   0.  +0.j        ,
1003              0.  +0.j        ,   0.  +0.j        ],
1004           [-12.5 +4.0614962j ,   0.  +0.j        ,   0.  +0.j        ,
1005              0.  +0.j        ,   0.  +0.j        ],
1006           [-12.5 -4.0614962j ,   0.  +0.j        ,   0.  +0.j        ,
1007              0.  +0.j        ,   0.  +0.j        ],
1008           [-12.5-17.20477401j,   0.  +0.j        ,   0.  +0.j        ,
1009              0.  +0.j        ,   0.  +0.j        ]])
1010
1011    """
1012    return _raw_fftnd(a, s, axes, fft, norm)
1013
1014
1015@array_function_dispatch(_fftn_dispatcher)
1016def ifft2(a, s=None, axes=(-2, -1), norm=None):
1017    """
1018    Compute the 2-dimensional inverse discrete Fourier Transform.
1019
1020    This function computes the inverse of the 2-dimensional discrete Fourier
1021    Transform over any number of axes in an M-dimensional array by means of
1022    the Fast Fourier Transform (FFT).  In other words, ``ifft2(fft2(a)) == a``
1023    to within numerical accuracy.  By default, the inverse transform is
1024    computed over the last two axes of the input array.
1025
1026    The input, analogously to `ifft`, should be ordered in the same way as is
1027    returned by `fft2`, i.e. it should have the term for zero frequency
1028    in the low-order corner of the two axes, the positive frequency terms in
1029    the first half of these axes, the term for the Nyquist frequency in the
1030    middle of the axes and the negative frequency terms in the second half of
1031    both axes, in order of decreasingly negative frequency.
1032
1033    Parameters
1034    ----------
1035    a : array_like
1036        Input array, can be complex.
1037    s : sequence of ints, optional
1038        Shape (length of each axis) of the output (``s[0]`` refers to axis 0,
1039        ``s[1]`` to axis 1, etc.).  This corresponds to `n` for ``ifft(x, n)``.
1040        Along each axis, if the given shape is smaller than that of the input,
1041        the input is cropped.  If it is larger, the input is padded with zeros.
1042        if `s` is not given, the shape of the input along the axes specified
1043        by `axes` is used.  See notes for issue on `ifft` zero padding.
1044    axes : sequence of ints, optional
1045        Axes over which to compute the FFT.  If not given, the last two
1046        axes are used.  A repeated index in `axes` means the transform over
1047        that axis is performed multiple times.  A one-element sequence means
1048        that a one-dimensional FFT is performed.
1049    norm : {"backward", "ortho", "forward"}, optional
1050        .. versionadded:: 1.10.0
1051
1052        Normalization mode (see `numpy.fft`). Default is "backward".
1053        Indicates which direction of the forward/backward pair of transforms
1054        is scaled and with what normalization factor.
1055
1056        .. versionadded:: 1.20.0
1057
1058            The "backward", "forward" values were added.
1059
1060    Returns
1061    -------
1062    out : complex ndarray
1063        The truncated or zero-padded input, transformed along the axes
1064        indicated by `axes`, or the last two axes if `axes` is not given.
1065
1066    Raises
1067    ------
1068    ValueError
1069        If `s` and `axes` have different length, or `axes` not given and
1070        ``len(s) != 2``.
1071    IndexError
1072        If an element of `axes` is larger than than the number of axes of `a`.
1073
1074    See Also
1075    --------
1076    numpy.fft : Overall view of discrete Fourier transforms, with definitions
1077         and conventions used.
1078    fft2 : The forward 2-dimensional FFT, of which `ifft2` is the inverse.
1079    ifftn : The inverse of the *n*-dimensional FFT.
1080    fft : The one-dimensional FFT.
1081    ifft : The one-dimensional inverse FFT.
1082
1083    Notes
1084    -----
1085    `ifft2` is just `ifftn` with a different default for `axes`.
1086
1087    See `ifftn` for details and a plotting example, and `numpy.fft` for
1088    definition and conventions used.
1089
1090    Zero-padding, analogously with `ifft`, is performed by appending zeros to
1091    the input along the specified dimension.  Although this is the common
1092    approach, it might lead to surprising results.  If another form of zero
1093    padding is desired, it must be performed before `ifft2` is called.
1094
1095    Examples
1096    --------
1097    >>> a = 4 * np.eye(4)
1098    >>> np.fft.ifft2(a)
1099    array([[1.+0.j,  0.+0.j,  0.+0.j,  0.+0.j], # may vary
1100           [0.+0.j,  0.+0.j,  0.+0.j,  1.+0.j],
1101           [0.+0.j,  0.+0.j,  1.+0.j,  0.+0.j],
1102           [0.+0.j,  1.+0.j,  0.+0.j,  0.+0.j]])
1103
1104    """
1105    return _raw_fftnd(a, s, axes, ifft, norm)
1106
1107
1108@array_function_dispatch(_fftn_dispatcher)
1109def rfftn(a, s=None, axes=None, norm=None):
1110    """
1111    Compute the N-dimensional discrete Fourier Transform for real input.
1112
1113    This function computes the N-dimensional discrete Fourier Transform over
1114    any number of axes in an M-dimensional real array by means of the Fast
1115    Fourier Transform (FFT).  By default, all axes are transformed, with the
1116    real transform performed over the last axis, while the remaining
1117    transforms are complex.
1118
1119    Parameters
1120    ----------
1121    a : array_like
1122        Input array, taken to be real.
1123    s : sequence of ints, optional
1124        Shape (length along each transformed axis) to use from the input.
1125        (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
1126        The final element of `s` corresponds to `n` for ``rfft(x, n)``, while
1127        for the remaining axes, it corresponds to `n` for ``fft(x, n)``.
1128        Along any axis, if the given shape is smaller than that of the input,
1129        the input is cropped.  If it is larger, the input is padded with zeros.
1130        if `s` is not given, the shape of the input along the axes specified
1131        by `axes` is used.
1132    axes : sequence of ints, optional
1133        Axes over which to compute the FFT.  If not given, the last ``len(s)``
1134        axes are used, or all axes if `s` is also not specified.
1135    norm : {"backward", "ortho", "forward"}, optional
1136        .. versionadded:: 1.10.0
1137
1138        Normalization mode (see `numpy.fft`). Default is "backward".
1139        Indicates which direction of the forward/backward pair of transforms
1140        is scaled and with what normalization factor.
1141
1142        .. versionadded:: 1.20.0
1143
1144            The "backward", "forward" values were added.
1145
1146    Returns
1147    -------
1148    out : complex ndarray
1149        The truncated or zero-padded input, transformed along the axes
1150        indicated by `axes`, or by a combination of `s` and `a`,
1151        as explained in the parameters section above.
1152        The length of the last axis transformed will be ``s[-1]//2+1``,
1153        while the remaining transformed axes will have lengths according to
1154        `s`, or unchanged from the input.
1155
1156    Raises
1157    ------
1158    ValueError
1159        If `s` and `axes` have different length.
1160    IndexError
1161        If an element of `axes` is larger than than the number of axes of `a`.
1162
1163    See Also
1164    --------
1165    irfftn : The inverse of `rfftn`, i.e. the inverse of the n-dimensional FFT
1166         of real input.
1167    fft : The one-dimensional FFT, with definitions and conventions used.
1168    rfft : The one-dimensional FFT of real input.
1169    fftn : The n-dimensional FFT.
1170    rfft2 : The two-dimensional FFT of real input.
1171
1172    Notes
1173    -----
1174    The transform for real input is performed over the last transformation
1175    axis, as by `rfft`, then the transform over the remaining axes is
1176    performed as by `fftn`.  The order of the output is as for `rfft` for the
1177    final transformation axis, and as for `fftn` for the remaining
1178    transformation axes.
1179
1180    See `fft` for details, definitions and conventions used.
1181
1182    Examples
1183    --------
1184    >>> a = np.ones((2, 2, 2))
1185    >>> np.fft.rfftn(a)
1186    array([[[8.+0.j,  0.+0.j], # may vary
1187            [0.+0.j,  0.+0.j]],
1188           [[0.+0.j,  0.+0.j],
1189            [0.+0.j,  0.+0.j]]])
1190
1191    >>> np.fft.rfftn(a, axes=(2, 0))
1192    array([[[4.+0.j,  0.+0.j], # may vary
1193            [4.+0.j,  0.+0.j]],
1194           [[0.+0.j,  0.+0.j],
1195            [0.+0.j,  0.+0.j]]])
1196
1197    """
1198    a = asarray(a)
1199    s, axes = _cook_nd_args(a, s, axes)
1200    a = rfft(a, s[-1], axes[-1], norm)
1201    for ii in range(len(axes)-1):
1202        a = fft(a, s[ii], axes[ii], norm)
1203    return a
1204
1205
1206@array_function_dispatch(_fftn_dispatcher)
1207def rfft2(a, s=None, axes=(-2, -1), norm=None):
1208    """
1209    Compute the 2-dimensional FFT of a real array.
1210
1211    Parameters
1212    ----------
1213    a : array
1214        Input array, taken to be real.
1215    s : sequence of ints, optional
1216        Shape of the FFT.
1217    axes : sequence of ints, optional
1218        Axes over which to compute the FFT.
1219    norm : {"backward", "ortho", "forward"}, optional
1220        .. versionadded:: 1.10.0
1221
1222        Normalization mode (see `numpy.fft`). Default is "backward".
1223        Indicates which direction of the forward/backward pair of transforms
1224        is scaled and with what normalization factor.
1225
1226        .. versionadded:: 1.20.0
1227
1228            The "backward", "forward" values were added.
1229
1230    Returns
1231    -------
1232    out : ndarray
1233        The result of the real 2-D FFT.
1234
1235    See Also
1236    --------
1237    rfftn : Compute the N-dimensional discrete Fourier Transform for real
1238            input.
1239
1240    Notes
1241    -----
1242    This is really just `rfftn` with different default behavior.
1243    For more details see `rfftn`.
1244
1245    """
1246    return rfftn(a, s, axes, norm)
1247
1248
1249@array_function_dispatch(_fftn_dispatcher)
1250def irfftn(a, s=None, axes=None, norm=None):
1251    """
1252    Computes the inverse of `rfftn`.
1253
1254    This function computes the inverse of the N-dimensional discrete
1255    Fourier Transform for real input over any number of axes in an
1256    M-dimensional array by means of the Fast Fourier Transform (FFT).  In
1257    other words, ``irfftn(rfftn(a), a.shape) == a`` to within numerical
1258    accuracy. (The ``a.shape`` is necessary like ``len(a)`` is for `irfft`,
1259    and for the same reason.)
1260
1261    The input should be ordered in the same way as is returned by `rfftn`,
1262    i.e. as for `irfft` for the final transformation axis, and as for `ifftn`
1263    along all the other axes.
1264
1265    Parameters
1266    ----------
1267    a : array_like
1268        Input array.
1269    s : sequence of ints, optional
1270        Shape (length of each transformed axis) of the output
1271        (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.). `s` is also the
1272        number of input points used along this axis, except for the last axis,
1273        where ``s[-1]//2+1`` points of the input are used.
1274        Along any axis, if the shape indicated by `s` is smaller than that of
1275        the input, the input is cropped.  If it is larger, the input is padded
1276        with zeros. If `s` is not given, the shape of the input along the axes
1277        specified by axes is used. Except for the last axis which is taken to
1278        be ``2*(m-1)`` where ``m`` is the length of the input along that axis.
1279    axes : sequence of ints, optional
1280        Axes over which to compute the inverse FFT. If not given, the last
1281        `len(s)` axes are used, or all axes if `s` is also not specified.
1282        Repeated indices in `axes` means that the inverse transform over that
1283        axis is performed multiple times.
1284    norm : {"backward", "ortho", "forward"}, optional
1285        .. versionadded:: 1.10.0
1286
1287        Normalization mode (see `numpy.fft`). Default is "backward".
1288        Indicates which direction of the forward/backward pair of transforms
1289        is scaled and with what normalization factor.
1290
1291        .. versionadded:: 1.20.0
1292
1293            The "backward", "forward" values were added.
1294
1295    Returns
1296    -------
1297    out : ndarray
1298        The truncated or zero-padded input, transformed along the axes
1299        indicated by `axes`, or by a combination of `s` or `a`,
1300        as explained in the parameters section above.
1301        The length of each transformed axis is as given by the corresponding
1302        element of `s`, or the length of the input in every axis except for the
1303        last one if `s` is not given.  In the final transformed axis the length
1304        of the output when `s` is not given is ``2*(m-1)`` where ``m`` is the
1305        length of the final transformed axis of the input.  To get an odd
1306        number of output points in the final axis, `s` must be specified.
1307
1308    Raises
1309    ------
1310    ValueError
1311        If `s` and `axes` have different length.
1312    IndexError
1313        If an element of `axes` is larger than than the number of axes of `a`.
1314
1315    See Also
1316    --------
1317    rfftn : The forward n-dimensional FFT of real input,
1318            of which `ifftn` is the inverse.
1319    fft : The one-dimensional FFT, with definitions and conventions used.
1320    irfft : The inverse of the one-dimensional FFT of real input.
1321    irfft2 : The inverse of the two-dimensional FFT of real input.
1322
1323    Notes
1324    -----
1325    See `fft` for definitions and conventions used.
1326
1327    See `rfft` for definitions and conventions used for real input.
1328
1329    The correct interpretation of the hermitian input depends on the shape of
1330    the original data, as given by `s`. This is because each input shape could
1331    correspond to either an odd or even length signal. By default, `irfftn`
1332    assumes an even output length which puts the last entry at the Nyquist
1333    frequency; aliasing with its symmetric counterpart. When performing the
1334    final complex to real transform, the last value is thus treated as purely
1335    real. To avoid losing information, the correct shape of the real input
1336    **must** be given.
1337
1338    Examples
1339    --------
1340    >>> a = np.zeros((3, 2, 2))
1341    >>> a[0, 0, 0] = 3 * 2 * 2
1342    >>> np.fft.irfftn(a)
1343    array([[[1.,  1.],
1344            [1.,  1.]],
1345           [[1.,  1.],
1346            [1.,  1.]],
1347           [[1.,  1.],
1348            [1.,  1.]]])
1349
1350    """
1351    a = asarray(a)
1352    s, axes = _cook_nd_args(a, s, axes, invreal=1)
1353    for ii in range(len(axes)-1):
1354        a = ifft(a, s[ii], axes[ii], norm)
1355    a = irfft(a, s[-1], axes[-1], norm)
1356    return a
1357
1358
1359@array_function_dispatch(_fftn_dispatcher)
1360def irfft2(a, s=None, axes=(-2, -1), norm=None):
1361    """
1362    Computes the inverse of `rfft2`.
1363
1364    Parameters
1365    ----------
1366    a : array_like
1367        The input array
1368    s : sequence of ints, optional
1369        Shape of the real output to the inverse FFT.
1370    axes : sequence of ints, optional
1371        The axes over which to compute the inverse fft.
1372        Default is the last two axes.
1373    norm : {"backward", "ortho", "forward"}, optional
1374        .. versionadded:: 1.10.0
1375
1376        Normalization mode (see `numpy.fft`). Default is "backward".
1377        Indicates which direction of the forward/backward pair of transforms
1378        is scaled and with what normalization factor.
1379
1380        .. versionadded:: 1.20.0
1381
1382            The "backward", "forward" values were added.
1383
1384    Returns
1385    -------
1386    out : ndarray
1387        The result of the inverse real 2-D FFT.
1388
1389    See Also
1390    --------
1391    rfft2 : The forward two-dimensional FFT of real input,
1392            of which `irfft2` is the inverse.
1393    rfft : The one-dimensional FFT for real input.
1394    irfft : The inverse of the one-dimensional FFT of real input.
1395    irfftn : Compute the inverse of the N-dimensional FFT of real input.
1396
1397    Notes
1398    -----
1399    This is really `irfftn` with different defaults.
1400    For more details see `irfftn`.
1401
1402    """
1403    return irfftn(a, s, axes, norm)
1404