This file is indexed.

/usr/include/relion-1.4/src/fftw.h is in librelion-dev-common 1.4+dfsg-1build1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
838
839
840
841
842
843
844
845
846
847
848
849
/***************************************************************************
 *
 * Author: "Sjors H.W. Scheres"
 * MRC Laboratory of Molecular Biology
 *
 * This program is free software; you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 2 of the License, or
 * (at your option) any later version.
 *
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * This complete copyright notice must be included in any revised version of the
 * source code. Additional authorship citations may be added, but existing
 * author citations must be preserved.
 ***************************************************************************/
/***************************************************************************
 *
 * Authors:    Roberto Marabini                 (roberto@cnb.csic.es)
 *             Carlos Oscar S. Sorzano          (coss@cnb.csic.es)
 *
 * Unidad de  Bioinformatica of Centro Nacional de Biotecnologia , CSIC
 *
 * This program is free software; you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 2 of the License, or
 * (at your option) any later version.
 *
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, write to the Free Software
 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
 * 02111-1307  USA
 *
 *  All comments concerning this program package may be sent to the
 *  e-mail address 'xmipp@cnb.csic.es'
 ***************************************************************************/

#ifndef __XmippFFTW_H
#define __XmippFFTW_H

#include <fftw3.h>
#include "src/multidim_array.h"
#include "src/funcs.h"
#include "src/tabfuncs.h"
#include "src/complex.h"

/** @defgroup FourierW FFTW Fourier transforms
  * @ingroup DataLibrary
  */

/** For all direct elements in the complex array in FFTW format.
 *
 * This macro is used to generate loops for the volume in an easy way. It
 * defines internal indexes 'k','i' and 'j' which ranges the volume using its
 * physical definition. It also defines 'kp', 'ip' and 'jp', which are the logical coordinates
 * It also works for 1D or 2D FFTW transforms
 *
 * @code
 * FOR_ALL_ELEMENTS_IN_FFTW_TRANSFORM(V)
 * {
 *     int r2 = jp*jp + ip*ip + kp*kp;
 *
 *     std::cout << "element at physical coords: "<< i<<" "<<j<<" "<<k<<" has value: "<<DIRECT_A3D_ELEM(m, k, i, j) << std::endl;
 *     std::cout << "its logical coords are: "<< ip<<" "<<jp<<" "<<kp<<std::endl;
 *     std::cout << "its distance from the origin = "<<sqrt(r2)<<std::endl;
 *
 * }
 * @endcode
 */
#define FOR_ALL_ELEMENTS_IN_FFTW_TRANSFORM(V) \
    for (long int k = 0, kp = 0; k<ZSIZE(V); k++, kp = (k < XSIZE(V)) ? k : k - ZSIZE(V)) \
    	for (long int i = 0, ip = 0 ; i<YSIZE(V); i++, ip = (i < XSIZE(V)) ? i : i - YSIZE(V)) \
    		for (long int j = 0, jp = 0; j<XSIZE(V); j++, jp = j)

/** For all direct elements in the complex array in FFTW format.
 *  The same as above, but now only for 2D images (this saves some time as k is not sampled
 */
#define FOR_ALL_ELEMENTS_IN_FFTW_TRANSFORM2D(V) \
	for (long int i = 0, ip = 0 ; i<YSIZE(V); i++, ip = (i < XSIZE(V)) ? i : i - YSIZE(V)) \
		for (long int j = 0, jp = 0; j<XSIZE(V); j++, jp = j)

/** FFTW volume element: Logical access.
 *
 * @code
 *
 * FFTW_ELEM(V, -1, -2, 1) = 1;
 * val = FFTW_ELEM(V, -1, -2, 1);
 * @endcode
 */
#define FFTW_ELEM(V, kp, ip, jp) \
    DIRECT_A3D_ELEM((V),((kp < 0) ? (kp + ZSIZE(V)) : (kp)), ((ip < 0) ? (ip + YSIZE(V)) : (ip)), (jp))

/** FFTW 2D image element: Logical access.
 *
 * @code
 *
 * FFTW2D_ELEM(V, --2, 1) = 1;
 * val = FFTW2D_ELEM(V, -2, 1);
 * @endcode
 */
#define FFTW2D_ELEM(V, ip, jp) \
    DIRECT_A2D_ELEM((V), ((ip < 0) ? (ip + YSIZE(V)) : (ip)), (jp))

/** Fourier Transformer class.
 * @ingroup FourierW
 *
 * The memory for the Fourier transform is handled by this object.
 * However, the memory for the real space image is handled externally
 * and this object only has a pointer to it.
 *
 * Here you have an example of use
 * @code
 * FourierTransformer transformer;
 * MultidimArray< Complex > Vfft;
 * transformer.FourierTransform(V(),Vfft,false);
 * MultidimArray<DOUBLE> Vmag;
 * Vmag.resize(Vfft);
 * FOR_ALL_ELEMENTS_IN_ARRAY3D(Vmag)
 *     Vmag(k,i,j)=20*log10(abs(Vfft(k,i,j)));
 * @endcode
 */
class FourierTransformer
{
public:
    /** Real array, in fact a pointer to the user array is stored. */
    MultidimArray<DOUBLE> *fReal;

     /** Complex array, in fact a pointer to the user array is stored. */
    MultidimArray<Complex > *fComplex;

    /** Fourier array  */
    MultidimArray< Complex > fFourier;

#ifdef FLOAT_PRECISION
    /* fftw Forward plan */
    fftwf_plan fPlanForward;

    /* fftw Backward plan */
    fftwf_plan fPlanBackward;
#else
    /* fftw Forward plan */
    fftw_plan fPlanForward;

    /* fftw Backward plan */
    fftw_plan fPlanBackward;
#endif

    /* number of threads*/
    int nthreads;

    /* Threads has been used in this program*/
    bool threadsSetOn;

// Public methods
public:
    /** Default constructor */
    FourierTransformer();

    /** Destructor */
    ~FourierTransformer();

    /** Copy constructor
     *
     * The created FourierTransformer is a perfect copy of the input array but with a
     * different memory assignment.
     *
     */
    FourierTransformer(const FourierTransformer& op);

    /** Set Number of threads
     * This function, which should be called once, performs any
     * one-time initialization required to use threads on your
     * system.
     *
     *  The nthreads argument indicates the number of threads you
     *  want FFTW to use (or actually, the maximum number). All
     *  plans subsequently created with any planner routine will use
     *  that many threads. You can call fftw_plan_with_nthreads,
     *  create some plans, call fftw_plan_with_nthreads again with a
     *  different argument, and create some more plans for a new
     *  number of threads. Plans already created before a call to
     *  fftw_plan_with_nthreads are unaffected. If you pass an
     *  nthreads argument of 1 (the default), threads are
     *  disabled for subsequent plans. */
    void setThreadsNumber(int tNumber);

    /** Compute the Fourier transform of a MultidimArray, 2D and 3D.
        If getCopy is false, an alias to the transformed data is returned.
        This is a faster option since a copy of all the data is avoided,
        but you need to be careful that an inverse Fourier transform may
        change the data.
        */
    template <typename T, typename T1>
        void FourierTransform(T& v, T1& V, bool getCopy=true)
        {
            setReal(v);
            Transform(FFTW_FORWARD);
            if (getCopy) getFourierCopy(V);
            else         getFourierAlias(V);
        }

    /** Compute the Fourier transform.
        The data is taken from the matrix with which the object was
        created. */
    void FourierTransform();

    /** Inforce Hermitian symmetry.
        If the Fourier transform risks of losing Hermitian symmetry,
        use this function to renforce it. */
    void enforceHermitianSymmetry();

    /** Compute the inverse Fourier transform.
        The result is stored in the same real data that was passed for
        the forward transform. The Fourier coefficients are taken from
        the internal Fourier coefficients */
    void inverseFourierTransform();

    /** Compute the inverse Fourier transform.
        New data is provided for the Fourier coefficients and the output
        can be any matrix1D, 2D or 3D. It is important that the output
        matrix is already resized to the right size before entering
        in this function. */
    template <typename T, typename T1>
        void inverseFourierTransform(T& V, T1& v)
        {
            setReal(v);
            setFourier(V);
            Transform(FFTW_BACKWARD);
        }

    /** Get Fourier coefficients. */
    template <typename T>
        void getFourierAlias(T& V) {V.alias(fFourier); return;}

    /** Get Fourier coefficients. */
    template <typename T>
        void getFourierCopy(T& V) {
            V.resize(fFourier);
            memcpy(MULTIDIM_ARRAY(V),MULTIDIM_ARRAY(fFourier),
                MULTIDIM_SIZE(fFourier)*2*sizeof(DOUBLE));
        }

    /** Return a complete Fourier transform (two halves).
    */
    template <typename T>
        void getCompleteFourier(T& V) {
            V.resize(*fReal);
            int ndim=3;
            if (ZSIZE(*fReal)==1)
            {
                ndim=2;
                if (YSIZE(*fReal)==1)
                    ndim=1;
            }
            switch (ndim)
            {
                case 1:
                    FOR_ALL_DIRECT_ELEMENTS_IN_ARRAY1D(V)
                        if (i<XSIZE(fFourier))
                            DIRECT_A1D_ELEM(V,i)=DIRECT_A1D_ELEM(fFourier,i);
                        else
                            DIRECT_A1D_ELEM(V,i)=
                                conj(DIRECT_A1D_ELEM(fFourier,
                                    XSIZE(*fReal)-i));
                    break;
                case 2:
                    FOR_ALL_DIRECT_ELEMENTS_IN_ARRAY2D(V)
                        if (j<XSIZE(fFourier))
                            DIRECT_A2D_ELEM(V,i,j)=
                                DIRECT_A2D_ELEM(fFourier,i,j);
                        else
                            DIRECT_A2D_ELEM(V,i,j)=
                                conj(DIRECT_A2D_ELEM(fFourier,
                                    (YSIZE(*fReal)-i)%YSIZE(*fReal),
                                     XSIZE(*fReal)-j));
                    break;
                case 3:
                    FOR_ALL_DIRECT_ELEMENTS_IN_ARRAY3D(V)
                        if (j<XSIZE(fFourier))
                            DIRECT_A3D_ELEM(V,k,i,j)=
                                DIRECT_A3D_ELEM(fFourier,k,i,j);
                        else
                            DIRECT_A3D_ELEM(V,k,i,j)=
                                conj(DIRECT_A3D_ELEM(fFourier,
                                    (ZSIZE(*fReal)-k)%ZSIZE(*fReal),
                                    (YSIZE(*fReal)-i)%YSIZE(*fReal),
                                     XSIZE(*fReal)-j));
                    break;
            }
        }

    /** Set one half of the FT in fFourier from the input complete Fourier transform (two halves).
        The fReal and fFourier already should have the right sizes
    */
    template <typename T>
        void setFromCompleteFourier(T& V) {
        int ndim=3;
        if (ZSIZE(*fReal)==1)
        {
            ndim=2;
            if (YSIZE(*fReal)==1)
                ndim=1;
        }
        switch (ndim)
        {
        case 1:
            FOR_ALL_DIRECT_ELEMENTS_IN_ARRAY1D(fFourier)
                DIRECT_A1D_ELEM(fFourier,i)=DIRECT_A1D_ELEM(V,i);
            break;
        case 2:
            FOR_ALL_DIRECT_ELEMENTS_IN_ARRAY2D(fFourier)
                DIRECT_A2D_ELEM(fFourier,i,j) = DIRECT_A2D_ELEM(V,i,j);
            break;
        case 3:
            FOR_ALL_DIRECT_ELEMENTS_IN_ARRAY3D(fFourier)
                DIRECT_A3D_ELEM(fFourier,k,i,j) = DIRECT_A3D_ELEM(V,k,i,j);
            break;
        }
    }

// Internal methods
public:
    /* Pointer to the array of DOUBLEs with which the plan was computed */
    DOUBLE * dataPtr;

    /* Pointer to the array of complex<DOUBLE> with which the plan was computed */
    Complex * complexDataPtr;

    /* Initialise all pointers to NULL */
    void init();

    /** Clear object */
    void clear();

    /** This calls fftw_cleanup.
     * NOTE!! When using multiple threads, only ONE thread can call this function, as it cleans up things that are shared among all threads...
 	 *  Therefore, this cleanup is something that needs to be done manually...
    */
    void cleanup();

    /** Destroy both forward and backward fftw planes (mutex locked */
    void destroyPlans();

    /** Computes the transform, specified in Init() function
        If normalization=true the forward transform is normalized
        (no normalization is made in the inverse transform)
        If normalize=false no normalization is performed and therefore
        the image is scaled by the number of pixels.
    */
    void Transform(int sign);

    /** Get the Multidimarray that is being used as input. */
    const MultidimArray<DOUBLE> &getReal() const;
    const MultidimArray<Complex > &getComplex() const;

    /** Set a Multidimarray for input.
        The data of img will be the one of fReal. In forward
        transforms it is not modified, but in backward transforms,
        the result will be stored in img. This means that the size
        of img cannot change between calls. */
    void setReal(MultidimArray<DOUBLE> &img);

    /** Set a Multidimarray for input.
        The data of img will be the one of fComplex. In forward
        transforms it is not modified, but in backward transforms,
        the result will be stored in img. This means that the size
        of img cannot change between calls. */
    void setReal(MultidimArray<Complex > &img);

    /** Set a Multidimarray for the Fourier transform.
        The values of the input array are copied in the internal array.
        It is assumed that the container for the real image as well as
        the one for the Fourier array are already resized.
        No plan is updated. */
    void setFourier(MultidimArray<Complex > &imgFourier);
};

// Randomize phases beyond the given shell (index)
void randomizePhasesBeyond(MultidimArray<DOUBLE> &I, int index);

/** Center an array, to have its origin at the origin of the FFTW
 *
 */
template <typename T>
void CenterFFT(MultidimArray< T >& v, bool forward)
{
    if ( v.getDim() == 1 )
    {
        // 1D
        MultidimArray< T > aux;
        int l, shift;

        l = XSIZE(v);
        aux.resize(l);
        shift = (int)(l / 2);

        if (!forward)
            shift = -shift;

        // Shift the input in an auxiliar vector
        for (int i = 0; i < l; i++)
        {
            int ip = i + shift;

            if (ip < 0)
                ip += l;
            else if (ip >= l)
                ip -= l;

            aux(ip) = DIRECT_A1D_ELEM(v, i);
        }

        // Copy the vector
        for (int i = 0; i < l; i++)
            DIRECT_A1D_ELEM(v, i) = DIRECT_A1D_ELEM(aux, i);
    }
    else if ( v.getDim() == 2 )
    {
        // 2D
        MultidimArray< T > aux;
        int l, shift;

        // Shift in the X direction
        l = XSIZE(v);
        aux.resize(l);
        shift = (int)(l / 2);

        if (!forward)
            shift = -shift;

        for (int i = 0; i < YSIZE(v); i++)
        {
            // Shift the input in an auxiliar vector
            for (int j = 0; j < l; j++)
            {
                int jp = j + shift;

                if (jp < 0)
                    jp += l;
                else if (jp >= l)
                    jp -= l;

                aux(jp) = DIRECT_A2D_ELEM(v, i, j);
            }

            // Copy the vector
            for (int j = 0; j < l; j++)
                DIRECT_A2D_ELEM(v, i, j) = DIRECT_A1D_ELEM(aux, j);
        }

        // Shift in the Y direction
        l = YSIZE(v);
        aux.resize(l);
        shift = (int)(l / 2);

        if (!forward)
            shift = -shift;

        for (int j = 0; j < XSIZE(v); j++)
        {
            // Shift the input in an auxiliar vector
            for (int i = 0; i < l; i++)
            {
                int ip = i + shift;

                if (ip < 0)
                    ip += l;
                else if (ip >= l)
                    ip -= l;

                aux(ip) = DIRECT_A2D_ELEM(v, i, j);
            }

            // Copy the vector
            for (int i = 0; i < l; i++)
                DIRECT_A2D_ELEM(v, i, j) = DIRECT_A1D_ELEM(aux, i);
        }
    }
    else if ( v.getDim() == 3 )
    {
        // 3D
        MultidimArray< T > aux;
        int l, shift;

        // Shift in the X direction
        l = XSIZE(v);
        aux.resize(l);
        shift = (int)(l / 2);

        if (!forward)
            shift = -shift;

        for (int k = 0; k < ZSIZE(v); k++)
            for (int i = 0; i < YSIZE(v); i++)
            {
                // Shift the input in an auxiliar vector
                for (int j = 0; j < l; j++)
                {
                    int jp = j + shift;

                    if (jp < 0)
                        jp += l;
                    else if (jp >= l)
                        jp -= l;

                    aux(jp) = DIRECT_A3D_ELEM(v, k, i, j);
                }

                // Copy the vector
                for (int j = 0; j < l; j++)
                    DIRECT_A3D_ELEM(v, k, i, j) = DIRECT_A1D_ELEM(aux, j);
            }

        // Shift in the Y direction
        l = YSIZE(v);
        aux.resize(l);
        shift = (int)(l / 2);

        if (!forward)
            shift = -shift;

        for (int k = 0; k < ZSIZE(v); k++)
            for (int j = 0; j < XSIZE(v); j++)
            {
                // Shift the input in an auxiliar vector
                for (int i = 0; i < l; i++)
                {
                    int ip = i + shift;

                    if (ip < 0)
                        ip += l;
                    else if (ip >= l)
                        ip -= l;

                    aux(ip) = DIRECT_A3D_ELEM(v, k, i, j);
                }

                // Copy the vector
                for (int i = 0; i < l; i++)
                    DIRECT_A3D_ELEM(v, k, i, j) = DIRECT_A1D_ELEM(aux, i);
            }

        // Shift in the Z direction
        l = ZSIZE(v);
        aux.resize(l);
        shift = (int)(l / 2);

        if (!forward)
            shift = -shift;

        for (int i = 0; i < YSIZE(v); i++)
            for (int j = 0; j < XSIZE(v); j++)
            {
                // Shift the input in an auxiliar vector
                for (int k = 0; k < l; k++)
                {
                    int kp = k + shift;
                    if (kp < 0)
                        kp += l;
                    else if (kp >= l)
                        kp -= l;

                    aux(kp) = DIRECT_A3D_ELEM(v, k, i, j);
                }

                // Copy the vector
                for (int k = 0; k < l; k++)
                    DIRECT_A3D_ELEM(v, k, i, j) = DIRECT_A1D_ELEM(aux, k);
            }
    }
    else
    {
    	v.printShape();
    	REPORT_ERROR("CenterFFT ERROR: Dimension should be 1, 2 or 3");
    }
}



// Window an FFTW-centered Fourier-transform to a given size
template<class T>
void windowFourierTransform(MultidimArray<T > &in,
			  			    MultidimArray<T > &out,
			  			    long int newdim)
{
	// Check size of the input array
	if (YSIZE(in) > 1 && YSIZE(in)/2 + 1 != XSIZE(in))
		REPORT_ERROR("windowFourierTransform ERROR: the Fourier transform should be of an image with equal sizes in all dimensions!");
	long int newhdim = newdim/2 + 1;

	// If same size, just return input
	if (newhdim == XSIZE(in))
	{
		out = in;
		return;
	}

	// Otherwise apply a windowing operation
	// Initialise output array
	switch (in.getDim())
	{
	case 1:
		out.initZeros(newhdim);
		break;
	case 2:
		out.initZeros(newdim, newhdim);
		break;
	case 3:
		out.initZeros(newdim, newdim, newhdim);
		break;
	default:
    	REPORT_ERROR("windowFourierTransform ERROR: dimension should be 1, 2 or 3!");
    }
	if (newhdim > XSIZE(in))
	{
		long int max_r2 = (XSIZE(in) -1) * (XSIZE(in) - 1);
		FOR_ALL_ELEMENTS_IN_FFTW_TRANSFORM(in)
		{
			// Make sure windowed FT has nothing in the corners, otherwise we end up with an asymmetric FT!
			if (kp*kp + ip*ip + jp*jp <= max_r2)
				FFTW_ELEM(out, kp, ip, jp) = FFTW_ELEM(in, kp, ip, jp);
		}
	}
	else
	{
		FOR_ALL_ELEMENTS_IN_FFTW_TRANSFORM(out)
		{
			FFTW_ELEM(out, kp, ip, jp) = FFTW_ELEM(in, kp, ip, jp);
		}
	}
}

// A resize operation in Fourier-space (i.e. changing the sampling of the Fourier Transform) by windowing in real-space
// If recenter=true, the real-space array will be recentered to have its origin at the origin of the FT
template<class T>
void resizeFourierTransform(MultidimArray<T > &in,
			  			    MultidimArray<T > &out,
			  			    long int newdim, bool do_recenter=true)
{
	// Check size of the input array
	if (YSIZE(in) > 1 && YSIZE(in)/2 + 1 != XSIZE(in))
		REPORT_ERROR("windowFourierTransform ERROR: the Fourier transform should be of an image with equal sizes in all dimensions!");
	long int newhdim = newdim/2 + 1;
	long int olddim = 2* (XSIZE(in) - 1);

	// If same size, just return input
	if (newhdim == XSIZE(in))
	{
		out = in;
		return;
	}

	// Otherwise apply a windowing operation
	MultidimArray<Complex > Fin;
	MultidimArray<DOUBLE> Min;
	FourierTransformer transformer;
	long int x0, y0, z0, xF, yF, zF;
	x0 = y0 = z0 = FIRST_XMIPP_INDEX(newdim);
	xF = yF = zF = LAST_XMIPP_INDEX(newdim);

	// Initialise output array
	switch (in.getDim())
	{
	case 1:
		Min.resize(olddim);
		y0=yF=z0=zF=0;
		break;
	case 2:
		Min.resize(olddim, olddim);
		z0=zF=0;
		break;
	case 3:
		Min.resize(olddim, olddim, olddim);
		break;
	default:
    	REPORT_ERROR("resizeFourierTransform ERROR: dimension should be 1, 2 or 3!");
    }

	// This is to handle DOUBLE-valued input arrays
	Fin.resize(ZSIZE(in), YSIZE(in), XSIZE(in));
	FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(in)
	{
		DIRECT_MULTIDIM_ELEM(Fin, n) = DIRECT_MULTIDIM_ELEM(in, n);
	}
	transformer.inverseFourierTransform(Fin, Min);
	Min.setXmippOrigin();
	if (do_recenter)
		CenterFFT(Min, false);

	// Now do the actual windowing in real-space
	Min.window(z0, y0, x0, zF, yF, xF);
	Min.setXmippOrigin();

	// If upsizing: mask the corners to prevent aliasing artefacts
	if (newdim > olddim)
	{
		FOR_ALL_ELEMENTS_IN_ARRAY3D(Min)
		{
			if (k*k + i*i + j*j > olddim*olddim/4)
			{
				A3D_ELEM(Min, k, i, j) = 0.;
			}
		}
	}

	// Recenter FFT back again
	if (do_recenter)
		CenterFFT(Min, true);

	// And do the inverse Fourier transform
	transformer.clear();
	transformer.FourierTransform(Min, out);
}

/** Fourier-Ring-Correlation between two multidimArrays using FFT
 * From precalculated Fourier Transforms
 * Simpler I/O than above.
 */
void getFSC(MultidimArray< Complex > &FT1,
		    MultidimArray< Complex > &FT2,
		    MultidimArray< DOUBLE > &fsc);

/** Fourier-Ring-Correlation between two multidimArrays using FFT
 * @ingroup FourierOperations
 * Simpler I/O than above.
 */
void getFSC(MultidimArray< DOUBLE > & m1,
		    MultidimArray< DOUBLE > & m2,
		    MultidimArray< DOUBLE > &fsc);

/** Scale matrix using Fourier transform
 * @ingroup FourierOperations
 * Ydim and Xdim define the output size, Mpmem is the matrix to scale
 */
//void selfScaleToSizeFourier(long int Ydim, long int Xdim, MultidimArray<DOUBLE>& Mpmem, int nthreads=1);

// Get precalculated AB-matrices for on-the-fly shift calculations
void getAbMatricesForShiftImageInFourierTransform(MultidimArray<Complex > &in,
									MultidimArray<Complex > &out,
									TabSine &tab_sin, TabCosine &tab_cos,
									DOUBLE oridim, DOUBLE shift_x, DOUBLE shift_y, DOUBLE shift_z = 0.);

// Shift an image through phase-shifts in its Fourier Transform
// Note that in and out may be the same array, in that case in is overwritten with the result
// if oridim is in pixels, xshift, yshift and zshift should be in pixels as well!
// or both can be in Angstroms
void shiftImageInFourierTransform(MultidimArray<Complex > &in,
								  MultidimArray<Complex > &out,
								  TabSine &tab_sin, TabCosine &tab_cos,
								  DOUBLE oridim, DOUBLE shift_x, DOUBLE shift_y, DOUBLE shift_z = 0.);

// Shift an image through phase-shifts in its Fourier Transform (without tabulated sine and cosine)
// Note that in and out may be the same array, in that case in is overwritten with the result
// if oridim is in pixels, xshift, yshift and zshift should be in pixels as well!
// or both can be in Angstroms
void shiftImageInFourierTransform(MultidimArray<Complex > &in,
						          MultidimArray<Complex > &out,
								  DOUBLE oridim, DOUBLE shift_x, DOUBLE shift_y, DOUBLE shift_z = 0.);

#define POWER_SPECTRUM 0
#define AMPLITUDE_SPECTRUM 1

/** Get the amplitude or power_class spectrum of the map in Fourier space.
 * @ingroup FourierOperations
    i.e. the radial average of the (squared) amplitudes of all Fourier components
*/
void getSpectrum(MultidimArray<DOUBLE> &Min,
                 MultidimArray<DOUBLE> &spectrum,
                 int spectrum_type=POWER_SPECTRUM);

/** Divide the input map in Fourier-space by the spectrum provided.
 * @ingroup FourierOperations
    If leave_origin_intact==true, the origin pixel will remain untouched
*/
void divideBySpectrum(MultidimArray<DOUBLE> &Min,
                      MultidimArray<DOUBLE> &spectrum,
                      bool leave_origin_intact=false);

/** Multiply the input map in Fourier-space by the spectrum provided.
 * @ingroup FourierOperations
    If leave_origin_intact==true, the origin pixel will remain untouched
*/
void multiplyBySpectrum(MultidimArray<DOUBLE> &Min,
                        MultidimArray<DOUBLE> &spectrum,
                        bool leave_origin_intact=false);

/** Perform a whitening of the amplitude/power_class spectrum of a 3D map
 * @ingroup FourierOperations
    If leave_origin_intact==true, the origin pixel will remain untouched
*/
void whitenSpectrum(MultidimArray<DOUBLE> &Min,
                    MultidimArray<DOUBLE> &Mout,
                    int spectrum_type=AMPLITUDE_SPECTRUM,
                    bool leave_origin_intact=false);

/** Adapts Min to have the same spectrum as spectrum_ref
 * @ingroup FourierOperations
    If only_amplitudes==true, the amplitude rather than the power_class spectrum will be equalized
*/
void adaptSpectrum(MultidimArray<DOUBLE> &Min,
                   MultidimArray<DOUBLE> &Mout,
                   const MultidimArray<DOUBLE> &spectrum_ref,
                   int spectrum_type=AMPLITUDE_SPECTRUM,
                   bool leave_origin_intact=false);

/** Kullback-Leibner divergence */
DOUBLE getKullbackLeibnerDivergence(MultidimArray<Complex > &Fimg,
		MultidimArray<Complex > &Fref, MultidimArray<DOUBLE> &sigma2,
		MultidimArray<DOUBLE> &p_i, MultidimArray<DOUBLE> &q_i,
		int highshell = -1, int lowshell = -1);


// Resize a map by windowing it's Fourier Transform
void resizeMap(MultidimArray<DOUBLE > &img, int newsize);

// Apply a B-factor to a map (given it's Fourier transform)
void applyBFactorToMap(MultidimArray<Complex > &FT, int ori_size, DOUBLE bfactor, DOUBLE angpix);

// Apply a B-factor to a map (given it's real-space array)
void applyBFactorToMap(MultidimArray<DOUBLE > &img, DOUBLE bfactor, DOUBLE angpix);

// Low-pass filter a map (given it's Fourier transform)
void lowPassFilterMap(MultidimArray<Complex > &FT, int ori_size,
		DOUBLE low_pass, DOUBLE angpix, int filter_edge_width = 2, bool do_highpass_instead = false);

// Low-pass and high-pass filter a map (given it's real-space array)
void lowPassFilterMap(MultidimArray<DOUBLE > &img, DOUBLE low_pass, DOUBLE angpix, int filter_edge_width = 2);
void highPassFilterMap(MultidimArray<DOUBLE > &img, DOUBLE low_pass, DOUBLE angpix, int filter_edge_width = 2);

/*
 *  Beamtilt x and y are given in mradians
 *  Wavelength in Angstrom, Cs in mm
 *  Phase shifts caused by the beamtilt will be calculated and applied to Fimg
 */
void selfApplyBeamTilt(MultidimArray<Complex > &Fimg, DOUBLE beamtilt_x, DOUBLE beamtilt_y,
		DOUBLE wavelength, DOUBLE Cs, DOUBLE angpix, int ori_size);

void applyBeamTilt(const MultidimArray<Complex > &Fin, MultidimArray<Complex > &Fout, DOUBLE beamtilt_x, DOUBLE beamtilt_y,
		DOUBLE wavelength, DOUBLE Cs, DOUBLE angpix, int ori_size);

#endif