File: xfftn.py

package info (click to toggle)
mpi4py-fft 2.0.6-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 720 kB
  • sloc: python: 3,053; ansic: 87; makefile: 42; sh: 33
file content (837 lines) | stat: -rw-r--r-- 27,897 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
#pylint: disable=no-name-in-module,unused-import
import numpy as np
from .factory import get_planned_FFT
from .utilities import FFTW_FORWARD, FFTW_BACKWARD, FFTW_REDFT00, FFTW_REDFT01, \
    FFTW_REDFT10, FFTW_REDFT11, FFTW_RODFT00, FFTW_RODFT01, FFTW_RODFT10, \
    FFTW_RODFT11, FFTW_MEASURE, FFTW_DESTROY_INPUT, FFTW_UNALIGNED, \
    FFTW_CONSERVE_MEMORY, FFTW_EXHAUSTIVE, FFTW_PRESERVE_INPUT, FFTW_PATIENT, \
    FFTW_ESTIMATE, FFTW_WISDOM_ONLY, C2C_FORWARD, C2C_BACKWARD, R2C, C2R, \
    FFTW_R2HC, FFTW_HC2R, FFTW_DHT, get_alignment, aligned, aligned_like

flag_dict = {key: val for key, val in locals().items()
             if key.startswith('FFTW_')}

dct_type = {
    1: FFTW_REDFT00,
    2: FFTW_REDFT10,
    3: FFTW_REDFT01,
    4: FFTW_REDFT11}

idct_type = {
    1: FFTW_REDFT00,
    2: FFTW_REDFT01,
    3: FFTW_REDFT10,
    4: FFTW_REDFT11}

dst_type = {
    1: FFTW_RODFT00,
    2: FFTW_RODFT10,
    3: FFTW_RODFT01,
    4: FFTW_RODFT11}

idst_type = {
    1: FFTW_RODFT00,
    2: FFTW_RODFT01,
    3: FFTW_RODFT10,
    4: FFTW_RODFT11}

def fftn(input_array, s=None, axes=(-1,), threads=1,
         flags=(FFTW_MEASURE,), output_array=None):
    """Return complex-to-complex forward transform object

    Parameters
    ----------
    input_array : complex array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the FFT.
    threads : int, optional
        Number of threads used in computing FFT.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : complex array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for complex-to-complex transforms

    Note
    ----
    This routine does not compute the fftn, it merely returns an instance of
    a class that can do it.
    The contents of the input_array may be overwritten during planning. Make
    sure to keep a copy if needed.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import fftn as plan_fftn
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='D')
    >>> fftn = plan_fftn(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = fftn()
    >>> print(B)
    [10.+0.j -2.+2.j -2.+0.j -2.-2.j]
    >>> assert id(A) == id(fftn.input_array)
    >>> assert id(B) == id(fftn.output_array)

    """
    kind = FFTW_FORWARD
    assert input_array.dtype.char in 'FDG'
    if output_array is None:
        n = get_alignment(input_array)
        output_array = aligned(input_array.shape, n,
                               input_array.dtype.char.upper())
    else:
        assert input_array.shape == output_array.shape
        assert output_array.dtype.char == input_array.dtype.char.upper()
    M = np.prod(np.take(input_array.shape, axes))
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, 1.0/M)

def ifftn(input_array, s=None, axes=(-1,), threads=1,
          flags=(FFTW_MEASURE,), output_array=None):
    """
    Return complex-to-complex inverse transform object

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the inverse FFT.
    threads : int, optional
        Number of threads used in computing FFT.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for complex-to-complex inverse
        transforms

    Note
    ----
    This routine does not compute the ifftn, it merely returns an instance of
    a class that can do it.
    The contents of the input_array may be overwritten during planning. Make
    sure that you keep a copy if needed.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import ifftn as plan_ifftn
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, FFTW_PRESERVE_INPUT, aligned
    >>> A = aligned(4, dtype='D')
    >>> ifftn = plan_ifftn(A, flags=(FFTW_ESTIMATE, FFTW_PRESERVE_INPUT))
    >>> A[:] = 1, 2, 3, 4
    >>> B = ifftn()
    >>> print(B)
    [10.+0.j -2.-2.j -2.+0.j -2.+2.j]
    >>> assert id(B) == id(ifftn.output_array)
    >>> assert id(A) == id(ifftn.input_array)

    """
    kind = FFTW_BACKWARD
    assert input_array.dtype.char in 'FDG'
    if output_array is None:
        output_array = aligned_like(input_array)
    else:
        assert input_array.shape == output_array.shape
    M = np.prod(np.take(input_array.shape, axes))
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, 1.0/M)

def rfftn(input_array, s=None, axes=(-1,), threads=1,
          flags=(FFTW_MEASURE,), output_array=None):
    """Return real-to-complex transform object

    Parameters
    ----------
    input_array : real array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the real to complex FFT.
    threads : int, optional
        Number of threads used in computing FFT.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for real-to-complex transforms

    Note
    ----
    This routine does not compute the rfftn, it merely returns an instance of
    a class that can do it.
    The contents of the input_array may be overwritten during planning. Make
    sure that you keep a copy if needed.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import rfftn as plan_rfftn
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='d')
    >>> rfftn = plan_rfftn(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = rfftn()
    >>> print(B)
    [10.+0.j -2.+2.j -2.+0.j]
    >>> assert id(A) == id(rfftn.input_array)
    >>> assert id(B) == id(rfftn.output_array)

    """
    kind = R2C
    assert input_array.dtype.char in 'fdg'
    if output_array is None:
        sz = list(input_array.shape)
        sz[axes[-1]] = input_array.shape[axes[-1]]//2+1
        dtype = input_array.dtype.char
        n = get_alignment(input_array)
        output_array = aligned(sz, n=n, dtype=np.dtype(dtype.upper()))
    else:
        assert input_array.shape[axes[-1]]//2+1 == output_array.shape[axes[-1]]
    M = np.prod(np.take(input_array.shape, axes))
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, 1.0/M)

def irfftn(input_array, s=None, axes=(-1,), threads=1,
           flags=(FFTW_MEASURE,), output_array=None):
    """Return inverse complex-to-real transform object

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Shape of output array along each of the transformed axes. Must be same
        length as axes (len(s) == len(axes)). If not given it is assumed that
        the shape of the output along the first transformed axis (i.e.,
        axes[-1]) is an even number. It is not possible to determine exactly,
        because for a real transform the output of a real array of length N is
        N//2+1. However, both N=4 and N=5 gives 4//2+1=3 and 5//2+1=3, so it is
        not possible to determine whether 4 or 5 is correct. Hence it must be
        given.
    axes : sequence of ints, optional
        Axes over which to compute the real to complex FFT.
    threads : int, optional
        Number of threads used in computing FFT.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for complex-to-real transforms

    Note
    ----
    This routine does not compute the irfftn, it merely returns an instance of
    a class that can do it.
    The irfftn is not possible to use with the FFTW_PRESERVE_INPUT flag.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import irfftn as plan_irfftn
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='D')
    >>> irfftn = plan_irfftn(A, flags=(FFTW_ESTIMATE,)) # no shape given for output
    >>> A[:] = 1, 2, 3, 4
    >>> B = irfftn()
    >>> print(B)
    [15. -4.  0. -1.  0. -4.]
    >>> irfftn = plan_irfftn(A, s=(7,), flags=(FFTW_ESTIMATE,)) # output shape given
    >>> B = irfftn()
    >>> print(B)
    [19.         -5.04891734 -0.30797853 -0.64310413 -0.64310413 -0.30797853
     -5.04891734]
    >>> assert id(B) == id(irfftn.output_array)
    >>> assert id(A) == id(irfftn.input_array)

    """
    kind = C2R
    assert input_array.dtype.char in 'FDG'
    assert FFTW_PRESERVE_INPUT not in flags
    sz = list(input_array.shape)
    if s is not None:
        assert len(axes) == len(s)
        for q, axis in zip(s, axes):
            sz[axis] = q
    else:
        sz[axes[-1]] = 2*sz[axes[-1]]-2
    if output_array is None:
        dtype = input_array.dtype.char
        n = get_alignment(input_array)
        output_array = aligned(sz, n=n, dtype=np.dtype(dtype.lower()))
    else:
        assert list(output_array.shape) == sz

    assert sz[axes[-1]]//2+1 == input_array.shape[axes[-1]]
    M = np.prod(np.take(output_array.shape, axes))
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, 1.0/M)

def dctn(input_array, s=None, axes=(-1,), type=2, threads=1,
         flags=(FFTW_MEASURE,), output_array=None):
    """Return discrete cosine transform object

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the real-to-real dct.
    type : int, optional
        Type of `dct <http://www.fftw.org/fftw3_doc/Real_002dto_002dReal-Transform-Kinds.html>`_

            - 1 - FFTW_REDFT00
            - 2 - FFTW_REDFT10,
            - 3 - FFTW_REDFT01,
            - 4 - FFTW_REDFT11
    threads : int, optional
        Number of threads used in computing dct.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for real-to-real dct transforms
        of given type

    Note
    ----
    This routine does not compute the dct, it merely returns an instance of
    a class that can do it.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import dctn as plan_dct
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='d')
    >>> dct = plan_dct(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = dct()
    >>> print(B)
    [20.         -6.30864406  0.         -0.44834153]
    >>> assert id(A) == id(dct.input_array)
    >>> assert id(B) == id(dct.output_array)

    """
    assert input_array.dtype.char in 'fdg'
    if output_array is None:
        output_array = aligned_like(input_array)
    else:
        assert input_array.shape == output_array.shape
    kind = dct_type[type]
    kind = [kind]*len(axes)
    M = get_normalization(kind, input_array.shape, axes)
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, M)

def idctn(input_array, s=None, axes=(-1,), type=2, threads=1,
          flags=(FFTW_MEASURE,), output_array=None):
    """Return inverse discrete cosine transform object

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the real-to-real idct.
    type : int, optional
        Type of `idct <http://www.fftw.org/fftw3_doc/Real_002dto_002dReal-Transform-Kinds.html>`_

            - 1 - FFTW_REDFT00
            - 2 - FFTW_REDFT01
            - 3 - FFTW_REDFT10
            - 4 - FFTW_REDFT11
    threads : int, optional
        Number of threads used in computing idct.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for real-to-real idct transforms
        of given type

    Note
    ----
    This routine does not compute the idct, it merely returns an instance of
    a class that can do it.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import idctn as plan_idct
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='d')
    >>> idct = plan_idct(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = idct()
    >>> print(B)
    [11.99962628 -9.10294322  2.61766184 -1.5143449 ]
    >>> assert id(A) == id(idct.input_array)
    >>> assert id(B) == id(idct.output_array)

    """
    assert input_array.dtype.char in 'fdg'
    if output_array is None:
        output_array = aligned_like(input_array)
    else:
        assert input_array.shape == output_array.shape
    kind = idct_type[type]
    kind = [kind]*len(axes)
    M = get_normalization(kind, input_array.shape, axes)
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, M)

def dstn(input_array, s=None, axes=(-1,), type=2, threads=1,
         flags=(FFTW_MEASURE,), output_array=None):
    """Return discrete sine transform object

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the real-to-real dst.
    type : int, optional
        Type of `dst <http://www.fftw.org/fftw3_doc/Real_002dto_002dReal-Transform-Kinds.html>`_

            - 1 - FFTW_RODFT00
            - 2 - FFTW_RODFT10
            - 3 - FFTW_RODFT01
            - 4 - FFTW_RODFT11
    threads : int, optional
        Number of threads used in computing dst.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for real-to-real dst transforms
        of given type

    Note
    ----
    This routine does not compute the dst, it merely returns an instance of
    a class that can do it.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import dstn as plan_dst
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='d')
    >>> dst = plan_dst(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = dst()
    >>> print(B)
    [13.06562965 -5.65685425  5.411961   -4.        ]
    >>> assert id(A) == id(dst.input_array)
    >>> assert id(B) == id(dst.output_array)

    """
    assert input_array.dtype.char in 'fdg'
    if output_array is None:
        output_array = aligned_like(input_array)
    else:
        assert input_array.shape == output_array.shape
    kind = dst_type[type]
    kind = [kind]*len(axes)
    M = get_normalization(kind, input_array.shape, axes)
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, M)

def idstn(input_array, s=None, axes=(-1,), type=2, threads=1,
          flags=(FFTW_MEASURE,), output_array=None):
    """Return inverse discrete sine transform object

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the real-to-real inverse dst.
    type : int, optional
        Type of `idst <http://www.fftw.org/fftw3_doc/Real_002dto_002dReal-Transform-Kinds.html>`_

            - 1 - FFTW_RODFT00
            - 2 - FFTW_RODFT01
            - 3 - FFTW_RODFT10
            - 4 - FFTW_RODFT11
    threads : int, optional
        Number of threads used in computing inverse dst.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for real-to-real idst transforms
        of given type

    Note
    ----
    This routine does not compute the idst, it merely returns an instance of
    a class that can do it.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import idstn as plan_idst
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='d')
    >>> idst = plan_idst(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = idst()
    >>> print(B)
    [13.13707118 -1.6199144   0.72323135 -0.51978306]
    >>> assert id(A) == id(idst.input_array)
    >>> assert id(B) == id(idst.output_array)

    """
    assert input_array.dtype.char in 'fdg'
    if output_array is None:
        output_array = aligned_like(input_array)
    else:
        assert input_array.shape == output_array.shape
    kind = idst_type[type]
    kind = [kind]*len(axes)
    M = get_normalization(kind, input_array.shape, axes)
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, M)

def ihfftn(input_array, s=None, axes=(-1,), threads=1,
           flags=(FFTW_MEASURE,), output_array=None):
    """Return inverse transform object for an array with Hermitian symmetry

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the ihfftn.
    threads : int, optional
        Number of threads used in computing ihfftn.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for real-to-complex ihfftn
        transforms

    Note
    ----
    This routine does not compute the ihfttn, it merely returns an instance of
    a class that can do it.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import ihfftn as plan_ihfftn
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='d')
    >>> ihfftn = plan_ihfftn(A, flags=(FFTW_ESTIMATE,))
    >>> A[:] = 1, 2, 3, 4
    >>> B = ihfftn()
    >>> print(B)
    [10.+0.j -2.+2.j -2.+0.j]
    >>> assert id(A) == id(ihfftn.input_array)
    >>> assert id(B) == id(ihfftn.output_array)

    """
    kind = R2C
    assert input_array.dtype.char in 'fdg'
    if output_array is None:
        dtype = input_array.dtype.char
        sz = list(input_array.shape)
        sz[axes[-1]] = input_array.shape[axes[-1]]//2+1
        n = get_alignment(input_array)
        output_array = aligned(sz, n=n, dtype=np.dtype(dtype.upper()))
    else:
        assert input_array.shape[axes[-1]]//2+1 == output_array.shape[axes[-1]]
    M = get_normalization(kind, input_array.shape, axes)
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, M)

def hfftn(input_array, s=None, axes=(-1,), threads=1,
          flags=(FFTW_MEASURE,), output_array=None):
    """Return transform object for an array with Hermitian symmetry

    Parameters
    ----------
    input_array : array
    s : sequence of ints, optional
        Not used - included for compatibility with Numpy
    axes : sequence of ints, optional
        Axes over which to compute the hfftn.
    threads : int, optional
        Number of threads used in computing hfftn.
    flags : sequence of ints, optional
        Flags from

            - FFTW_MEASURE
            - FFTW_EXHAUSTIVE
            - FFTW_PATIENT
            - FFTW_DESTROY_INPUT
            - FFTW_PRESERVE_INPUT
            - FFTW_UNALIGNED
            - FFTW_CONSERVE_MEMORY
            - FFTW_ESTIMATE
    output_array : array, optional
        Array to be used as output array. Must be of correct shape, type,
        strides and alignment

    Returns
    -------
    :class:`.fftwf_xfftn.FFT`, :class:`.fftw_xfftn.FFT` or :class:`.fftwl_xfftn.FFT`
        An instance of the return type configured for complex-to-real hfftn
        transforms

    Note
    ----
    This routine does not compute the hfttn, it merely returns an instance of
    a class that can do it.

    Examples
    --------
    >>> import numpy as np
    >>> from mpi4py_fft.fftw import hfftn as plan_hfftn
    >>> from mpi4py_fft.fftw import FFTW_ESTIMATE, aligned
    >>> A = aligned(4, dtype='D')
    >>> hfftn = plan_hfftn(A, flags=(FFTW_ESTIMATE,)) # no shape given for output
    >>> A[:] = 1, 2, 3, 4
    >>> B = hfftn()
    >>> print(B)
    [15. -4.  0. -1.  0. -4.]
    >>> hfftn = plan_hfftn(A, s=(7,), flags=(FFTW_ESTIMATE,)) # output shape given
    >>> B = hfftn()
    >>> print(B)
    [19.         -5.04891734 -0.30797853 -0.64310413 -0.64310413 -0.30797853
     -5.04891734]
    >>> assert id(B) == id(hfftn.output_array)
    >>> assert id(A) == id(hfftn.input_array)

    """
    kind = C2R
    assert input_array.dtype.char in 'FDG'
    sz = list(input_array.shape)
    if s is not None:
        assert len(axes) == len(s)
        for q, axis in zip(s, axes):
            sz[axis] = q
    else:
        sz[axes[-1]] = 2*sz[axes[-1]]-2
    if output_array is None:
        dtype = input_array.dtype.char
        n = get_alignment(input_array)
        output_array = aligned(sz, n=n, dtype=np.dtype(dtype.lower()))
    else:
        assert list(output_array.shape) == sz
    assert sz[axes[-1]]//2+1 == input_array.shape[axes[-1]]
    M = get_normalization(kind, sz, axes)
    return get_planned_FFT(input_array, output_array, axes, kind, threads,
                           flags, M)

def get_normalization(kind, shape, axes):
    """Return normalization factor for multidimensional transform

    The normalization factor is, for Fourier transforms::

        1./np.prod(np.take(shape, axes))

    where shape is the global shape of the array that is input to the
    forward transform, and axes are the axes transformed over.

    For real-to-real transforms the normalization factor for each axis is

        - REDFT00 - 2(N-1)
        - REDFT01 - 2N
        - REDFT10 - 2N
        - REDFT11 - 2N
        - RODFT00 - 2(N+1)
        - RODFT01 - 2N
        - RODFT10 - 2N
        - RODFT11 - 2N

    where N is the length of the input array along that axis.

    Parameters
    ----------
    kind : sequence of ints
        The kind of transform along each axis
    shape : sequence of ints
        The shape of the global transformed array (input to the forward
        transform)
    axes : sequence of ints
        The axes transformed over

    Note
    ----
    The returned normalization factor is the *inverse* of the product of the
    normalization factors for the axes it is transformed over.

    """
    kind = [kind]*len(axes) if isinstance(kind, int) else kind
    assert len(kind) == len(axes)
    M = 1
    for knd, axis in zip(kind, axes):
        N = shape[axis]
        if knd == FFTW_RODFT00:
            M *= 2*(N+1)
        elif knd == FFTW_REDFT00:
            M *= 2*(N-1)
        elif knd in (FFTW_RODFT01, FFTW_RODFT10, FFTW_RODFT11,
                     FFTW_REDFT01, FFTW_REDFT10, FFTW_REDFT11):
            M *= 2*N
        else:
            M *= N
    return 1./M

inverse = {
    FFTW_RODFT11: FFTW_RODFT11,
    FFTW_REDFT11: FFTW_REDFT11,
    FFTW_RODFT01: FFTW_RODFT10,
    FFTW_RODFT10: FFTW_RODFT01,
    FFTW_REDFT01: FFTW_REDFT10,
    FFTW_REDFT10: FFTW_REDFT01,
    FFTW_RODFT00: FFTW_RODFT00,
    FFTW_REDFT00: FFTW_REDFT00,
    rfftn: irfftn,
    irfftn: rfftn,
    fftn: ifftn,
    ifftn: fftn,
    dctn: idctn,
    idctn: dctn,
    dstn: idstn,
    idstn: dstn,
    hfftn: ihfftn,
    ihfftn: hfftn
}