GammaLib 2.0.0
Loading...
Searching...
No Matches
GFft.cpp
Go to the documentation of this file.
1/***************************************************************************
2 * GFft.hpp - Fast Fourier transformation class *
3 * ----------------------------------------------------------------------- *
4 * copyright (C) 2016-2020 by Juergen Knoedlseder *
5 * ----------------------------------------------------------------------- *
6 * *
7 * This program is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * This program is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have received a copy of the GNU General Public License *
18 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
19 * *
20 ***************************************************************************/
21/**
22 * @file GFft.cpp
23 * @brief Fast Fourier transformation class implementation
24 * @author Juergen Knoedlseder
25 */
26
27/* __ Includes ___________________________________________________________ */
28#ifdef HAVE_CONFIG_H
29#include <config.h>
30#endif
31#include "GTools.hpp"
32#include "GMath.hpp"
33#include "GException.hpp"
34#include "GFft.hpp"
35#include "GNdarray.hpp"
36
37/* __ Method name definitions ____________________________________________ */
38#define G_OP_ADD "GFft::operator+=(GFft&)"
39#define G_OP_SUB "GFft::operator-=(GFft&)"
40#define G_OP_MUL "GFft::operator*=(GFft&)"
41#define G_OP_DIV "GFft::operator/=(GFft&)"
42#define G_FORWARD "GFft::forward(GNdarray&)"
43#define G_BACKWARD "GFft::backward()"
44
45
46/*==========================================================================
47 = =
48 = Constructors/destructors =
49 = =
50 ==========================================================================*/
51
52/***********************************************************************//**
53 * @brief Void constructor
54 ***************************************************************************/
56{
57 // Initialise class members
59
60 // Return
61 return;
62}
63
64
65/***********************************************************************//**
66 * @brief N-dimensional array constructor
67 *
68 * @param[in] array N-dimensional array.
69 *
70 * Constructs a Fast Fourier Transformation from a n-dimensional array.
71 ***************************************************************************/
72GFft::GFft(const GNdarray& array)
73{
74 // Initialise class members
76
77 // Perform foward transformation
78 forward(array);
79
80 // Return
81 return;
82}
83
84
85/***********************************************************************//**
86 * @brief Copy constructor
87 *
88 * @param[in] fft Fast Fourier Transform.
89 ***************************************************************************/
90GFft::GFft(const GFft& fft)
91{
92 // Initialise class members
94
95 // Copy members
96 copy_members(fft);
97
98 // Return
99 return;
100}
101
102
103/***********************************************************************//**
104 * @brief Destructor
105 ***************************************************************************/
107{
108 // Free members
109 free_members();
110
111 // Return
112 return;
113}
114
115
116/*==========================================================================
117 = =
118 = Operators =
119 = =
120 ==========================================================================*/
121
122/***********************************************************************//**
123 * @brief Assignment operator
124 *
125 * @param[in] fft Fast Fourier Transform.
126 * @return Fast Fourier Transform.
127 ***************************************************************************/
129{
130 // Execute only if object is not identical
131 if (this != &fft) {
132
133 // Free members
134 free_members();
135
136 // Initialise private members
137 init_members();
138
139 // Copy members
140 copy_members(fft);
141
142 } // endif: object was not identical
143
144 // Return this object
145 return *this;
146}
147
148
149/***********************************************************************//**
150 * @brief Unary addition operator
151 *
152 * @param[in] fft Fast Fourier Transform.
153 * @return Fast Fourier Transform.
154 *
155 * Adds a Fast Fourier Transform to another Fast Fourier Transform. The
156 * method requires that both FFTs have the same dimension.
157 ***************************************************************************/
159{
160 // Throw an exception if the FFTs have not the same shape
162
163 // Add elements
164 for (int i = 0; i < m_size; ++i) {
165 *(m_data+i) += *(fft.m_data+i);
166 }
167
168 // Return Fast Fourier Transform
169 return *this;
170}
171
172
173/***********************************************************************//**
174 * @brief Unary subtraction operator
175 *
176 * @param[in] fft Fast Fourier transform.
177 * @return Fast Fourier transform.
178 *
179 * Subtracts a Fast Fourier Transform from another Fast Fourier Transform.
180 * The method requires that both FFTs have the same dimension.
181 ***************************************************************************/
183{
184 // Throw an exception if the FFTs have not the same shape
186
187 // Subtract elements
188 for (int i = 0; i < m_size; ++i) {
189 *(m_data+i) -= *(fft.m_data+i);
190 }
191
192 // Return Fast Fourier Transform
193 return *this;
194}
195
196
197/***********************************************************************//**
198 * @brief Unary multiplication operator
199 *
200 * @param[in] fft Fast Fourier transform.
201 * @return Fast Fourier transform.
202 *
203 * Multiplies a Fast Fourier Transform to another Fast Fourier Transform.
204 * The method requires that both FFTs have the same dimension.
205 ***************************************************************************/
207{
208 // Throw an exception if the FFTs have not the same shape
210
211 // Multiply elements
212 for (int i = 0; i < m_size; ++i) {
213 *(m_data+i) *= *(fft.m_data+i);
214 }
215
216 // Return Fast Fourier Transform
217 return *this;
218}
219
220
221/***********************************************************************//**
222 * @brief Unary division operator
223 *
224 * @param[in] fft Fast Fourier transform.
225 * @return Fast Fourier transform.
226 *
227 * Divides a Fast Fourier Transform by another Fast Fourier Transform.
228 * The method requires that both FFTs have the same dimension.
229 ***************************************************************************/
231{
232 // Throw an exception if the FFTs have not the same shape
234
235 // Divide elements
236 for (int i = 0; i < m_size; ++i) {
237 *(m_data+i) /= *(fft.m_data+i);
238 }
239
240 // Return Fast Fourier Transform
241 return *this;
242}
243
244
245/***********************************************************************//**
246 * @brief Unary minus operator
247 *
248 * @return Fast Fourier transformation.
249 *
250 * Negates all elements of the Fast Fourier Transformation.
251 ***************************************************************************/
253{
254 // Copy FFT
255 GFft fft = *this;
256
257 // Negate all elements
258 for (int i = 0; i < fft.size(); ++i) {
259 fft(i) = -fft(i);
260 }
261
262 // Return FFT
263 return fft;
264}
265
266
267/*==========================================================================
268 = =
269 = Public methods =
270 = =
271 ==========================================================================*/
272
273/***********************************************************************//**
274 * @brief Clear Fast Fourier Transform
275 ***************************************************************************/
276void GFft::clear(void)
277{
278 // Free members
279 free_members();
280
281 // Initialise private members
282 init_members();
283
284 // Return
285 return;
286}
287
288
289/***********************************************************************//**
290 * @brief Clone Fast Fourier Transform
291 *
292 * @return Pointer to deep copy of Fast Fourier Transform.
293 ***************************************************************************/
294GFft* GFft::clone(void) const
295{
296 // Clone FFT
297 return new GFft(*this);
298}
299
300
301/***********************************************************************//**
302 * @brief Forward Fast Fourier Transform
303 *
304 * @param[in] array N-dimensional array.
305 *
306 * GException::feature_not_implemented
307 * Array dimension is not supported.
308 *
309 * Performs discret Fast Fourier forward transformation for 1- and
310 * 2-dimensional arrays.
311 *
312 * For one dimensional arrays of length \f$N\f$ the transformation is given
313 * by
314 *
315 * \f[
316 * X_k = \sum_{n=0}^{N-1} z_n e^{-i 2 \pi \frac{k n}{N}}
317 * \f]
318 *
319 * and for two dimensional arrays of length \f$ M \times N\f$ the
320 * transformation is given by
321 *
322 * \f[
323 * X_{k,l} = \sum_{n=0}^{N-1} \left(
324 * \sum_{m=0}^{M-1} z_{m,n} e^{-i 2 \pi \frac{m k}{M}} \right)
325 * e^{-i 2 \pi \frac{n l}{N}}
326 * \f]
327 ***************************************************************************/
328void GFft::forward(const GNdarray& array)
329{
330 // Set data
331 set_data(array);
332
333 // Continue only if there are elements in the data member
334 if (m_size > 0) {
335
336 // Case A: array is 1-dimensional
337 if (array.dim() == 1) {
338 transform(m_data, 1, m_shape[0], m_wavetable[0], true);
339 }
340
341 // Case B: array is 2-dimensional
342 else if (array.dim() == 2) {
343
344 // Loop over rows
345 for (int row = 0; row < m_shape[1]; ++row) {
346 transform(m_data + row*m_strides[1], m_strides[0], m_shape[0],
347 m_wavetable[0], true);
348 }
349
350 // Loop over columns
351 for (int col = 0; col < m_shape[0]; ++col) {
352 transform(m_data + col*m_strides[0], m_strides[1], m_shape[1],
353 m_wavetable[1], true);
354 }
355
356 }
357
358 // ... otherwise throw an exception
359 else {
360 std::string msg = "Fast Fourier Transformation for "+
361 gammalib::str(array.dim())+"-dimensional array "
362 "not implemented.";
364 }
365
366 } // endif: there were elements in the FFT data member
367
368 // Return
369 return;
370}
371
372
373/***********************************************************************//**
374 * @brief Backward Fast Fourier Transform
375 *
376 * @return N-dimensional array.
377 *
378 * GException::feature_not_implemented
379 * Array dimension is not supported.
380 *
381 * Performs the inverse (or backward) discret Fast Fourier forward
382 * transformation for 1- and 2-dimensional arrays.
383 *
384 * For one dimensional arrays of length \f$N\f$ the transformation is given
385 * by
386 *
387 * \f[
388 * z_k = \frac{1}{N} \sum_{n=0}^{N-1} X_n e^{i 2 \pi \frac{k n}{N}}
389 * \f]
390 *
391 * and for two dimensional arrays of length \f$ M \times N\f$ the
392 * transformation is given by
393 *
394 * \f[
395 * z_{k,l} = \frac{1}{N} \sum_{n=0}^{N-1} \left(
396 * \frac{1}{M} \sum_{m=0}^{M-1} X_{m,n}
397 * e^{i 2 \pi \frac{m k}{M}} \right)
398 * e^{i 2 \pi \frac{n l}{N}}
399 * \f]
400 ***************************************************************************/
402{
403 // Initialise n-dimensional array
404 GNdarray array(shape());
405
406 // Continue only if there are elements in the array
407 if (m_size > 0) {
408
409 // Create copy of complex array
410 std::complex<double>* data = new std::complex<double>[m_size];
411 for (int i = 0; i < m_size; ++i) {
412 *(data+i) = *(m_data+i);
413 }
414
415 // Case A: array is 1-dimensional
416 if (array.dim() == 1) {
417
418 // Perform backward transformation (circumvent const correctness)
419 const_cast<GFft*>(this)->transform(data, 1, m_shape[0], m_wavetable[0], false);
420
421 // Get normalisation factor
422 const double norm = 1.0 / double(m_shape[0]);
423
424 // Extract 1-dimensional array and normalise with 1/n
425 for (int i = 0; i < m_size; ++i) {
426 array(i) = (data+i)->real() * norm;
427 }
428
429 }
430
431 // Case B: array is 2-dimensional
432 else if (array.dim() == 2) {
433
434 // Loop over rows
435 for (int row = 0; row < m_shape[1]; ++row) {
436 const_cast<GFft*>(this)->transform(data + row*m_strides[1], m_strides[0], m_shape[0],
437 m_wavetable[0], false);
438 }
439
440 // Loop over columns
441 for (int col = 0; col < m_shape[0]; ++col) {
442 const_cast<GFft*>(this)->transform(data + col*m_strides[0], m_strides[1], m_shape[1],
443 m_wavetable[1], false);
444 }
445
446 // Get normalisation factor
447 const double norm = 1.0 / double(m_shape[0]*m_shape[1]);
448
449 // Extract 2-dimensional array and normalise with 1/(n*m)
450 for (int i = 0; i < m_size; ++i) {
451 array(i) = (data+i)->real() * norm;
452 }
453
454 }
455
456 // ... otherwise throw an exception
457 else {
458
459 // Delete copy of complex array
460 delete [] data;
461
462 // Throw exception
463 std::string msg = "Fast Fourier Transformation for "+
464 gammalib::str(array.dim())+"-dimensional array "
465 "not implemented.";
467 }
468
469 // Delete copy of complex array
470 delete [] data;
471
472 } // endif: there were elements in the array
473
474 // Return N-dimensional array
475 return array;
476}
477
478
479/***********************************************************************//**
480 * @brief Print Fast Fourier Transform information
481 *
482 * @param[in] chatter Chattiness.
483 * @return String containing Fast Fourier Transform information.
484 ***************************************************************************/
485std::string GFft::print(const GChatter& chatter) const
486{
487 // Initialise result string
488 std::string result;
489
490 // Continue only if chatter is not silent
491 if (chatter != SILENT) {
492
493 // Append header
494 result.append("=== GFft ===");
495
496 // Append array information
497 result.append("\n"+gammalib::parformat("Dimension"));
498 result.append(gammalib::str(dim()));
499 result.append("\n"+gammalib::parformat("Shape"));
500 result.append("(");
501 for (int i = 0; i < dim(); ++i) {
502 if (i > 0) {
503 result.append(",");
504 }
505 result.append(gammalib::str(m_shape[i]));
506 }
507 result.append(")");
508 result.append("\n"+gammalib::parformat("Size"));
509 result.append(gammalib::str(m_size));
510
511 // VERBOSE: Put all FFT elements in string
512 if (chatter >= VERBOSE) {
513 result.append("\n"+gammalib::parformat("Elements"));
514 for (int i = 0; i < size(); ++i) {
515 if (i > 0) {
516 result.append(",");
517 }
518 result.append(gammalib::str((*this)(i)));
519 }
520 result.append(")");
521 }
522
523 } // endif: chatter was not silent
524
525 // Return result
526 return result;
527}
528
529
530/*==========================================================================
531 = =
532 = Private methods =
533 = =
534 ==========================================================================*/
535
536/***********************************************************************//**
537 * @brief Initialise class members
538 ***************************************************************************/
540{
541 // Initialise members
542 m_size = 0;
543 m_data = NULL;
544 m_shape.clear();
545 m_strides.clear();
546 m_wavetable.clear();
547
548 // Return
549 return;
550}
551
552
553/***********************************************************************//**
554 * @brief Copy class members
555 *
556 * @param[in] fft Fast Fourier Transform.
557 ***************************************************************************/
558void GFft::copy_members(const GFft& fft)
559{
560 // Copy members
561 m_size = fft.m_size;
562 m_shape = fft.m_shape;
563 m_strides = fft.m_strides;
565
566 // Copy data
567 if (m_size > 0) {
568 m_data = new std::complex<double>[m_size];
569 for (int i = 0; i < m_size; ++i) {
570 *(m_data+i) = *(fft.m_data+i);
571 }
572 }
573
574 // Return
575 return;
576}
577
578
579/***********************************************************************//**
580 * @brief Delete class members
581 ***************************************************************************/
583{
584 // Delete data
585 if (m_data != NULL) {
586 delete [] m_data;
587 }
588
589 // Reset size and pointer
590 m_size = 0;
591 m_data = NULL;
592
593 // Return
594 return;
595}
596
597
598/***********************************************************************//**
599 * @brief Set data from n-dimensional array
600 *
601 * @param[in] array N-dimensional array
602 ***************************************************************************/
603void GFft::set_data(const GNdarray& array)
604{
605 // Clear information
606 free_members();
607 init_members();
608
609 // Continue only if there are elements in the array
610 if (array.size() > 0) {
611
612 // Set size of array
613 m_size = array.size();
614
615 // Allocate data and initialise data to zero
616 m_data = new std::complex<double>[m_size];
617
618 // Set real part of data from n-dimensional array
619 for (int i = 0; i < m_size; ++i) {
620 *(m_data+i) = std::complex<double>(array(i), 0.0);
621 }
622
623 // Save shape and strides
624 m_shape = array.shape();
625 m_strides = array.strides();
626
627 // Set trigonometric coefficients
628 for (int i = 0; i < m_shape.size(); ++i) {
629 GFftWavetable wavetable(m_shape[i]);
630 m_wavetable.push_back(wavetable);
631 }
632
633 } // endif: there were elements in the array
634
635 // Return
636 return;
637}
638
639
640/***********************************************************************//**
641 * @brief Check if FFT has the same shape
642 *
643 * @param[in] fft Fast Fourier Transform.
644 * @return True if the FFT has the same shape.
645 ***************************************************************************/
646bool GFft::has_same_shape(const GFft& fft) const
647{
648 // Check if the array has the same dimension
649 bool identity = m_shape.size() == fft.m_shape.size();
650
651 // Check if the shape is identical. Break as soon as one shape value
652 // differs.
653 if (identity) {
654 for (int i = 0; i < m_shape.size(); ++i) {
655 if (m_shape[i] != fft.m_shape[i]) {
656 identity = false;
657 break;
658 }
659 }
660 }
661
662 // Return result
663 return identity;
664}
665
666
667/***********************************************************************//**
668 * @brief Throw exception if FFT shapes differ
669 *
670 * @param[in] method Method that throws exception.
671 * @param[in] fft Fast Fourier Transform.
672 *
673 * @exception GException::invalid_argument
674 * FFT shapes differ.
675 ***************************************************************************/
676void GFft::require_same_shape(const std::string& method,
677 const GFft& fft) const
678{
679 // If the shape differs then throw an exception
680 if (!has_same_shape(fft)) {
681
682 // Compose exception message
683 std::string msg = "Incompatible FFT dimensions (";
684 for (int i = 0; i < m_shape.size(); ++i) {
685 if (i > 0) {
686 msg += ", ";
687 }
688 msg += gammalib::str(m_shape[i]);
689 }
690 msg += ") <=> (";
691 for (int i = 0; i < fft.m_shape.size(); ++i) {
692 if (i > 0) {
693 msg += ", ";
694 }
695 msg += gammalib::str(fft.m_shape[i]);
696 }
697 msg += ").";
698
699 // Throw exception
700 throw GException::invalid_argument(method, msg);
701
702 } // endif: array shapes differed
703
704 // Return
705 return;
706}
707
708
709/***********************************************************************//**
710 * @brief Perform Fast Fourier Transform
711 *
712 * @param[in] data Pointer to complex array to be transformed.
713 * @param[in] stride Step size when traversing complex array.
714 * @param[in] n Number of elements in complex array.
715 * @param[in] wavetable Trigonometric coefficients.
716 * @param[in] forward Forward transform if true, otherwise backward transform.
717 ***************************************************************************/
718void GFft::transform(std::complex<double>* data,
719 const int& stride,
720 const int& n,
721 const GFftWavetable& wavetable,
722 const bool& forward)
723{
724 // Allocate scratch space and initialize to zero
725 std::complex<double>* scratch = new std::complex<double>[n];
726 std::complex<double> zero(0.0, 0.0);
727 for (int i = 0; i < n; ++i) {
728 *(scratch+i) = zero;
729 }
730
731 // Initialise transformation state
732 int state = 0;
733
734 // Initialise factorisation product
735 int product = 1;
736
737 // Set sign
738 int sign = (forward) ? -1 : +1;
739
740 // Loop over all factors
741 for (int i = 0; i < wavetable.factors(); ++i) {
742
743 // Get factorisation factor and start index
744 int factor = wavetable.factor(i);
745 int index = wavetable.index(i);
746
747 // Update factorisation product
748 product *= factor;
749
750 // Set state dependent stuff
751 std::complex<double>* in;
752 std::complex<double>* out;
753 int istride;
754 int ostride;
755 if (state == 0) {
756 in = data;
757 istride = stride;
758 out = scratch;
759 ostride = 1;
760 state = 1;
761 }
762 else {
763 in = scratch;
764 istride = 1;
765 out = data;
766 ostride = stride;
767 state = 0;
768 }
769
770 // Call factor dependent method
771 if (factor == 2) {
772 factor2(in, istride, out, ostride, wavetable, sign, product, n, index);
773 }
774 else if (factor == 3) {
775 factor3(in, istride, out, ostride, wavetable, sign, product, n, index);
776 }
777 else if (factor == 4) {
778 factor4(in, istride, out, ostride, wavetable, sign, product, n, index);
779 }
780 else if (factor == 5) {
781 factor5(in, istride, out, ostride, wavetable, sign, product, n, index);
782 }
783 else if (factor == 6) {
784 factor6(in, istride, out, ostride, wavetable, sign, product, n, index);
785 }
786 else if (factor == 7) {
787 factor7(in, istride, out, ostride, wavetable, sign, product, n, index);
788 }
789 else {
790 factorn(in, istride, out, ostride, wavetable, sign, factor, product, n, index);
791 }
792
793 } // endfor: looped over all factors
794
795 // In case the loop was exited in state 1 the results are in the scratch
796 // array and we need to copy the results back from the scratch array to
797 // the data array
798 if (state == 1) {
799 for (int i = 0; i < n; ++i) {
800 *(data+stride*i) = *(scratch+i);
801 }
802 }
803
804 // Free scratch space
805 delete [] scratch;
806
807 // Return
808 return;
809}
810
811
812/***********************************************************************//**
813 * @brief Compute FFT for factor 2
814 *
815 * @param[in] in Pointer to input array.
816 * @param[in] istride Step size when traversing input array.
817 * @param[in] out Pointer to output array.
818 * @param[in] ostride Step size when traversing output array.
819 * @param[in] wavetable Trigonometric coefficients.
820 * @param[in] sign Forward (-1) or backward (+1).
821 * @param[in] product ...
822 * @param[in] n Logical array length.
823 * @param[in] index Start index of trigonometric coefficients.
824 ***************************************************************************/
825void GFft::factor2(const std::complex<double>* in,
826 const int& istride,
827 std::complex<double>* out,
828 const int& ostride,
829 const GFftWavetable& wavetable,
830 const int& sign,
831 const int& product,
832 const int& n,
833 const int& index)
834{
835 // Compute ...
836 const int factor = 2;
837 const int m = n / factor;
838 const int q = n / product;
839 const int p_1 = product / factor;
840 const int jump = (factor - 1) * p_1;
841
842 // Loop over ...
843 for (int k = 0, i = 0, j = 0; k < q; ++k, j += jump) {
844
845 // Extract coefficients from wavetable
846 std::vector<std::complex<double> > w = get_w(wavetable, index, k, q, 1, sign);
847
848 // Compute x = W(2) z
849 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
850
851 // Get z
852 const std::complex<double> z0 = *(in+istride*i);
853 const std::complex<double> z1 = *(in+istride*(i+m));
854
855 // Compute x
856 const std::complex<double> x0 = z0 + z1;
857 const std::complex<double> x1 = z0 - z1;
858
859 // out = w * x
860 *(out+ostride*j) = x0;
861 *(out+ostride*(j+p_1)) = w[0] * x1;
862
863 } // endfor: k1
864
865 } // endfor: k
866
867 // Return
868 return;
869}
870
871
872/***********************************************************************//**
873 * @brief Compute FFT for factor 3
874 *
875 * @param[in] in Pointer to input array.
876 * @param[in] istride Step size when traversing input array.
877 * @param[in] out Pointer to output array.
878 * @param[in] ostride Step size when traversing output array.
879 * @param[in] wavetable Trigonometric coefficients.
880 * @param[in] sign Forward (-1) or backward (+1).
881 * @param[in] product ...
882 * @param[in] n Logical array length.
883 * @param[in] index Start index of trigonometric coefficients.
884 ***************************************************************************/
885void GFft::factor3(const std::complex<double>* in,
886 const int& istride,
887 std::complex<double>* out,
888 const int& ostride,
889 const GFftWavetable& wavetable,
890 const int& sign,
891 const int& product,
892 const int& n,
893 const int& index)
894{
895 // Compute ...
896 const int factor = 3;
897 const int m = n / factor;
898 const int q = n / product;
899 const int p_1 = product / factor;
900 const int jump = (factor - 1) * p_1;
901
902
903 // Precompute some factors
904 const double tau = std::sqrt(3.0) / 2.0;
905
906 // Loop over ...
907 for (int k = 0, i = 0, j = 0; k < q; ++k, j += jump) {
908
909 // Extract coefficients from wavetable
910 std::vector<std::complex<double> > w = get_w(wavetable, index, k, q, 2, sign);
911
912 // Compute x = W(3) z
913 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
914
915 // Get z
916 const std::complex<double> z0 = *(in+istride*i);
917 const std::complex<double> z1 = *(in+istride*(i+m));
918 const std::complex<double> z2 = *(in+istride*(i+2*m));
919
920 // Compute t
921 const std::complex<double> t1 = z1 + z2;
922 const std::complex<double> t2 = z0 - t1/2.0;
923 const std::complex<double> t3 = double((int)sign) * tau * (z1 - z2);
924
925 // Compute x
926 const std::complex<double> x0 = z0 + t1;
927 const std::complex<double> x1 = t2 + timesi(t3);
928 const std::complex<double> x2 = t2 - timesi(t3);
929
930 // out = w * x
931 *(out+ostride*j) = x0;
932 *(out+ostride*(j+p_1)) = w[0] * x1;
933 *(out+ostride*(j+2*p_1)) = w[1] * x2;
934
935 } // endfor: k1
936
937 } // endfor: k
938
939 // Return
940 return;
941}
942
943
944/***********************************************************************//**
945 * @brief Compute FFT for factor 4
946 *
947 * @param[in] in Pointer to input array.
948 * @param[in] istride Step size when traversing input array.
949 * @param[in] out Pointer to output array.
950 * @param[in] ostride Step size when traversing output array.
951 * @param[in] wavetable Trigonometric coefficients.
952 * @param[in] sign Forward (-1) or backward (+1).
953 * @param[in] product ...
954 * @param[in] n Logical array length.
955 * @param[in] index Start index of trigonometric coefficients.
956 ***************************************************************************/
957void GFft::factor4(const std::complex<double>* in,
958 const int& istride,
959 std::complex<double>* out,
960 const int& ostride,
961 const GFftWavetable& wavetable,
962 const int& sign,
963 const int& product,
964 const int& n,
965 const int& index)
966{
967 // Compute ...
968 const int factor = 4;
969 const int m = n / factor;
970 const int q = n / product;
971 const int p_1 = product / factor;
972 const int jump = (factor - 1) * p_1;
973
974 // Loop over ...
975 for (int k = 0, i = 0, j = 0; k < q; ++k, j += jump) {
976
977 // Extract coefficients from wavetable
978 std::vector<std::complex<double> > w = get_w(wavetable, index, k, q, 3, sign);
979
980 // Compute x = W(4) z
981 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
982
983 // Get z
984 const std::complex<double> z0 = *(in+istride*i);
985 const std::complex<double> z1 = *(in+istride*(i+m));
986 const std::complex<double> z2 = *(in+istride*(i+2*m));
987 const std::complex<double> z3 = *(in+istride*(i+3*m));
988
989 // Compute t
990 const std::complex<double> t1 = z0 + z2;
991 const std::complex<double> t2 = z1 + z3;
992 const std::complex<double> t3 = z0 - z2;
993 const std::complex<double> t4 = double((int)sign) * (z1 - z3);
994
995 // Compute x
996 const std::complex<double> x0 = t1 + t2;
997 const std::complex<double> x1 = t3 + timesi(t4);
998 const std::complex<double> x2 = t1 - t2;
999 const std::complex<double> x3 = t3 - timesi(t4);
1000
1001 // out = w * x
1002 *(out+ostride*j) = x0;
1003 *(out+ostride*(j+p_1)) = w[0] * x1;
1004 *(out+ostride*(j+2*p_1)) = w[1] * x2;
1005 *(out+ostride*(j+3*p_1)) = w[2] * x3;
1006
1007 } // endfor: k1
1008
1009 } // endfor: k
1010
1011 // Return
1012 return;
1013}
1014
1015
1016/***********************************************************************//**
1017 * @brief Compute FFT for factor 5
1018 *
1019 * @param[in] in Pointer to input array.
1020 * @param[in] istride Step size when traversing input array.
1021 * @param[in] out Pointer to output array.
1022 * @param[in] ostride Step size when traversing output array.
1023 * @param[in] wavetable Trigonometric coefficients.
1024 * @param[in] sign Forward (-1) or backward (+1).
1025 * @param[in] product ...
1026 * @param[in] n Logical array length.
1027 * @param[in] index Start index of trigonometric coefficients.
1028 ***************************************************************************/
1029void GFft::factor5(const std::complex<double>* in,
1030 const int& istride,
1031 std::complex<double>* out,
1032 const int& ostride,
1033 const GFftWavetable& wavetable,
1034 const int& sign,
1035 const int& product,
1036 const int& n,
1037 const int& index)
1038{
1039 // Compute ...
1040 const int factor = 5;
1041 const int m = n / factor;
1042 const int q = n / product;
1043 const int p_1 = product / factor;
1044 const int jump = (factor - 1) * p_1;
1045
1046 // Precompute some factors
1047 const double sin_2pi_by_5 = std::sin(gammalib::twopi / 5.0);
1048 const double sin_2pi_by_10 = std::sin(gammalib::twopi / 10.0);
1049 const double tau = std::sqrt(5.0) / 4.0;
1050
1051 // Loop over ...
1052 for (int k = 0, i = 0, j = 0; k < q; ++k, j += jump) {
1053
1054 // Extract coefficients from wavetable
1055 std::vector<std::complex<double> > w = get_w(wavetable, index, k, q, 4, sign);
1056
1057 // Compute x = W(5) z
1058 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
1059
1060 // Get z
1061 const std::complex<double> z0 = *(in+istride*i);
1062 const std::complex<double> z1 = *(in+istride*(i+m));
1063 const std::complex<double> z2 = *(in+istride*(i+2*m));
1064 const std::complex<double> z3 = *(in+istride*(i+3*m));
1065 const std::complex<double> z4 = *(in+istride*(i+4*m));
1066
1067 // Compute t
1068 const std::complex<double> t1 = z1 + z4;
1069 const std::complex<double> t2 = z2 + z3;
1070 const std::complex<double> t3 = z1 - z4;
1071 const std::complex<double> t4 = z2 - z3;
1072 const std::complex<double> t5 = t1 + t2;
1073 const std::complex<double> t6 = tau * (t1 - t2);
1074 const std::complex<double> t7 = z0 - t5/4.0;
1075 const std::complex<double> t8 = t7 + t6;
1076 const std::complex<double> t9 = t7 - t6;
1077 const std::complex<double> t10 = double((int)sign) *
1078 (sin_2pi_by_5*t3 + sin_2pi_by_10*t4);
1079 const std::complex<double> t11 = double((int)sign) *
1080 (sin_2pi_by_10*t3 - sin_2pi_by_5*t4);
1081
1082 // Compute x
1083 const std::complex<double> x0 = z0 + t5;
1084 const std::complex<double> x1 = t8 + timesi(t10);
1085 const std::complex<double> x2 = t9 + timesi(t11);
1086 const std::complex<double> x3 = t9 - timesi(t11);
1087 const std::complex<double> x4 = t8 - timesi(t10);
1088
1089 // Compute out = w * x
1090 *(out+ostride*j) = x0;
1091 *(out+ostride*(j+p_1)) = w[0] * x1;
1092 *(out+ostride*(j+2*p_1)) = w[1] * x2;
1093 *(out+ostride*(j+3*p_1)) = w[2] * x3;
1094 *(out+ostride*(j+4*p_1)) = w[3] * x4;
1095
1096 } // endfor: k1
1097
1098 } // endfor: k
1099
1100 // Return
1101 return;
1102}
1103
1104
1105/***********************************************************************//**
1106 * @brief Compute FFT for factor 6
1107 *
1108 * @param[in] in Pointer to input array.
1109 * @param[in] istride Step size when traversing input array.
1110 * @param[in] out Pointer to output array.
1111 * @param[in] ostride Step size when traversing output array.
1112 * @param[in] wavetable Trigonometric coefficients.
1113 * @param[in] sign Forward (-1) or backward (+1).
1114 * @param[in] product ...
1115 * @param[in] n Logical array length.
1116 * @param[in] index Start index of trigonometric coefficients.
1117 ***************************************************************************/
1118void GFft::factor6(const std::complex<double>* in,
1119 const int& istride,
1120 std::complex<double>* out,
1121 const int& ostride,
1122 const GFftWavetable& wavetable,
1123 const int& sign,
1124 const int& product,
1125 const int& n,
1126 const int& index)
1127{
1128 // Compute ...
1129 const int factor = 6;
1130 const int m = n / factor;
1131 const int q = n / product;
1132 const int p_1 = product / factor;
1133 const int jump = (factor - 1) * p_1;
1134
1135 // Precompute some factors
1136 const double tau = std::sqrt(3.0) / 2.0;
1137
1138 // Loop over ...
1139 for (int k = 0, i = 0, j = 0; k < q; ++k, j += jump) {
1140
1141 // Extract coefficients from wavetable
1142 std::vector<std::complex<double> > w = get_w(wavetable, index, k, q, 5, sign);
1143
1144 // Compute x = W(6) z. W(6) is a combination of sums and differences of
1145 // W(3) acting on the even and odd elements of z
1146 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
1147
1148 // Get z
1149 const std::complex<double> z0 = *(in+istride*i);
1150 const std::complex<double> z1 = *(in+istride*(i+m));
1151 const std::complex<double> z2 = *(in+istride*(i+2*m));
1152 const std::complex<double> z3 = *(in+istride*(i+3*m));
1153 const std::complex<double> z4 = *(in+istride*(i+4*m));
1154 const std::complex<double> z5 = *(in+istride*(i+5*m));
1155
1156 // Compute ta
1157 const std::complex<double> ta1 = z2 + z4;
1158 const std::complex<double> ta2 = z0 - ta1/2.0;
1159 const std::complex<double> ta3 = double((int)sign) * tau * (z2 - z4);
1160
1161 // Compute a
1162 const std::complex<double> a0 = z0 + ta1;
1163 const std::complex<double> a1 = ta2 + timesi(ta3);
1164 const std::complex<double> a2 = ta2 - timesi(ta3);
1165
1166 // Compute tb
1167 const std::complex<double> tb1 = z5 + z1;
1168 const std::complex<double> tb2 = z3 - tb1/2.0;
1169 const std::complex<double> tb3 = double((int)sign) * tau * (z5 - z1);
1170
1171 // Compute b
1172 const std::complex<double> b0 = z3 + tb1;
1173 const std::complex<double> b1 = tb2 + timesi(tb3);
1174 const std::complex<double> b2 = tb2 - timesi(tb3);
1175
1176 // Compute x
1177 const std::complex<double> x0 = a0 + b0;
1178 const std::complex<double> x1 = a1 - b1;
1179 const std::complex<double> x2 = a2 + b2;
1180 const std::complex<double> x3 = a0 - b0;
1181 const std::complex<double> x4 = a1 + b1;
1182 const std::complex<double> x5 = a2 - b2;
1183
1184 // Compute out = w * x
1185 *(out+ostride*j) = x0;
1186 *(out+ostride*(j+p_1)) = w[0] * x1;
1187 *(out+ostride*(j+2*p_1)) = w[1] * x2;
1188 *(out+ostride*(j+3*p_1)) = w[2] * x3;
1189 *(out+ostride*(j+4*p_1)) = w[3] * x4;
1190 *(out+ostride*(j+5*p_1)) = w[4] * x5;
1191
1192 } // endfor: k1
1193
1194 } // endfor: k
1195
1196 // Return
1197 return;
1198}
1199
1200
1201/***********************************************************************//**
1202 * @brief Compute FFT for factor 7
1203 *
1204 * @param[in] in Pointer to input array.
1205 * @param[in] istride Step size when traversing input array.
1206 * @param[in] out Pointer to output array.
1207 * @param[in] ostride Step size when traversing output array.
1208 * @param[in] wavetable Trigonometric coefficients.
1209 * @param[in] sign Forward (-1) or backward (+1).
1210 * @param[in] product ...
1211 * @param[in] n Logical array length.
1212 * @param[in] index Start index of trigonometric coefficients.
1213 ***************************************************************************/
1214void GFft::factor7(const std::complex<double>* in,
1215 const int& istride,
1216 std::complex<double>* out,
1217 const int& ostride,
1218 const GFftWavetable& wavetable,
1219 const int& sign,
1220 const int& product,
1221 const int& n,
1222 const int& index)
1223{
1224 // Compute ...
1225 const int factor = 7;
1226 const int m = n / factor;
1227 const int q = n / product;
1228 const int p_1 = product / factor;
1229 const int jump = (factor - 1) * p_1;
1230
1231 // Precompute some factors
1232 static const double twopi7 = gammalib::twopi / 7.0;
1233 static const double c1 = std::cos(1.0 * twopi7);
1234 static const double c2 = std::cos(2.0 * twopi7);
1235 static const double c3 = std::cos(3.0 * twopi7);
1236 static const double s1 = std::sin(1.0 * twopi7);
1237 static const double s2 = std::sin(2.0 * twopi7);
1238 static const double s3 = std::sin(3.0 * twopi7);
1239 static const double tau1 = (c1 + c2 + c3) / 3.0 - 1.0;
1240 static const double tau2 = (2.0 * c1 - c2 - c3) / 3.0;
1241 static const double tau3 = (c1 - 2.0*c2 + c3) / 3.0;
1242 static const double tau4 = (c1 + c2 - 2.0 * c3) / 3.0;
1243 static const double tau5 = (s1 + s2 - s3) / 3.0;
1244 static const double tau6 = (2.0 * s1 - s2 + s3) / 3.0;
1245 static const double tau7 = (s1 - 2.0 * s2 - s3) / 3.0;
1246 static const double tau8 = (s1 + s2 + 2.0 * s3) / 3.0;
1247
1248 // Loop over ...
1249 for (int k = 0, i = 0, j = 0; k < q; ++k, j += jump) {
1250
1251 // Extract coefficients from wavetable
1252 std::vector<std::complex<double> > w = get_w(wavetable, index, k, q, 6, sign);
1253
1254 // Compute x = W(7) z
1255 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
1256
1257 // Get z
1258 const std::complex<double> z0 = *(in+istride*i);
1259 const std::complex<double> z1 = *(in+istride*(i+m));
1260 const std::complex<double> z2 = *(in+istride*(i+2*m));
1261 const std::complex<double> z3 = *(in+istride*(i+3*m));
1262 const std::complex<double> z4 = *(in+istride*(i+4*m));
1263 const std::complex<double> z5 = *(in+istride*(i+5*m));
1264 const std::complex<double> z6 = *(in+istride*(i+6*m));
1265
1266 // Compute t
1267 const std::complex<double> t0 = z1 + z6;
1268 const std::complex<double> t1 = z1 - z6;
1269 const std::complex<double> t2 = z2 + z5;
1270 const std::complex<double> t3 = z2 - z5;
1271 const std::complex<double> t4 = z4 + z3;
1272 const std::complex<double> t5 = z4 - z3;
1273 const std::complex<double> t6 = t2 + t0;
1274 const std::complex<double> t7 = t5 + t3;
1275
1276 // Compute b
1277 const std::complex<double> b0 = z0 + t6 + t4;
1278 const std::complex<double> b1 = tau1 * (t6 + t4);
1279 const std::complex<double> b2 = tau2 * (t0 - t4);
1280 const std::complex<double> b3 = tau3 * (t4 - t2);
1281 const std::complex<double> b4 = tau4 * (t2 - t0);
1282 const std::complex<double> b5 = double(-(int)sign) * tau5 * (t7 + t1);
1283 const std::complex<double> b6 = double(-(int)sign) * tau6 * (t1 - t5);
1284 const std::complex<double> b7 = double(-(int)sign) * tau7 * (t5 - t3);
1285 const std::complex<double> b8 = double(-(int)sign) * tau8 * (t3 - t1);
1286
1287 // Compute T
1288 const std::complex<double> T0 = b0 + b1;
1289 const std::complex<double> T1 = b2 + b3;
1290 const std::complex<double> T2 = b4 - b3;
1291 const std::complex<double> T3 = -b2 - b4;
1292 const std::complex<double> T4 = b6 + b7;
1293 const std::complex<double> T5 = b8 - b7;
1294 const std::complex<double> T6 = -b8 - b6;
1295 const std::complex<double> T7 = T0 + T1;
1296 const std::complex<double> T8 = T0 + T2;
1297 const std::complex<double> T9 = T0 + T3;
1298 const std::complex<double> T10 = T4 + b5;
1299 const std::complex<double> T11 = T5 + b5;
1300 const std::complex<double> T12 = T6 + b5;
1301
1302 // Compute x
1303 const std::complex<double> x0 = b0;
1304 const std::complex<double> x1 = T7 - timesi(T10);
1305 const std::complex<double> x2 = T9 - timesi(T12);
1306 const std::complex<double> x3 = T8 + timesi(T11);
1307 const std::complex<double> x4 = T8 - timesi(T11);
1308 const std::complex<double> x5 = T9 + timesi(T12);
1309 const std::complex<double> x6 = T7 + timesi(T10);
1310
1311 // Compute out = w * x
1312 *(out+ostride*j) = x0;
1313 *(out+ostride*(j+p_1)) = w[0] * x1;
1314 *(out+ostride*(j+2*p_1)) = w[1] * x2;
1315 *(out+ostride*(j+3*p_1)) = w[2] * x3;
1316 *(out+ostride*(j+4*p_1)) = w[3] * x4;
1317 *(out+ostride*(j+5*p_1)) = w[4] * x5;
1318 *(out+ostride*(j+6*p_1)) = w[5] * x6;
1319
1320 } // endfor: k1
1321
1322 } // endfor: k
1323
1324 // Return
1325 return;
1326}
1327
1328
1329/***********************************************************************//**
1330 * @brief Compute FFT for arbitrary factor
1331 *
1332 * @param[in] in Pointer to input array.
1333 * @param[in] istride Step size when traversing input array.
1334 * @param[in] out Pointer to output array.
1335 * @param[in] ostride Step size when traversing output array.
1336 * @param[in] wavetable Trigonometric coefficients.
1337 * @param[in] sign Forward (-1) or backward (+1).
1338 * @param[in] factor Factor.
1339 * @param[in] product ...
1340 * @param[in] n Logical array length.
1341 * @param[in] index Start index of trigonometric coefficients.
1342 ***************************************************************************/
1343void GFft::factorn(std::complex<double>* in,
1344 const int& istride,
1345 std::complex<double>* out,
1346 const int& ostride,
1347 const GFftWavetable& wavetable,
1348 const int& sign,
1349 const int& factor,
1350 const int& product,
1351 const int& n,
1352 const int& index)
1353{
1354 // Compute ...
1355 const int m = n / factor;
1356 const int q = n / product;
1357 const int p_1 = product / factor;
1358 const int jump = (factor - 1) * p_1;
1359
1360 // ...
1361 for (int i = 0; i < m; ++i) {
1362 *(out+ostride*i) = *(in+istride*i);
1363 }
1364
1365 // ...
1366 for (int e = 1; e < (factor - 1) / 2 + 1; ++e) {
1367 for (int i = 0; i < m; ++i) {
1368 const int idx = i + e * m;
1369 const int idxc = i + (factor - e) * m;
1370 *(out+ostride*idx) = *(in+istride*idx) + *(in+istride*idxc);
1371 *(out+ostride*idxc) = *(in+istride*idx) - *(in+istride*idxc);
1372 }
1373 }
1374
1375 // e = 0
1376 for (int i = 0; i < m; ++i) {
1377 *(in+istride*i) = *(out+ostride*i);
1378 }
1379 for (int e = 1; e < (factor - 1) / 2 + 1; ++e) {
1380 for (int i = 0; i < m; ++i) {
1381 const int idx = i + e * m;
1382 *(in+istride*i) += *(out+ostride*idx);
1383 }
1384 }
1385
1386 // ...
1387 for (int e = 1; e < (factor-1)/2 + 1; ++e) {
1388
1389 //
1390 int idx = e * q;
1391 const int idx_step = e * q;
1392
1393 //
1394 double w_real, w_imag;
1395
1396 //
1397 const int em = e * m;
1398 const int ecm = (factor - e) * m;
1399
1400 // ...
1401 for (int i = 0; i < m; ++i) {
1402 *(in+istride*(i+em)) = *(out+ostride*i);
1403 *(in+istride*(i+ecm)) = *(out+ostride*i);
1404 }
1405
1406 // ...
1407 for (int e1 = 1; e1 < (factor - 1) / 2 + 1; ++e1) {
1408
1409 //
1410 if (idx == 0) {
1411 w_real = 1.0;
1412 w_imag = 0.0;
1413 }
1414 else {
1415
1416 // Compute indices
1417 int twiddle = index + idx - 1;
1418
1419 // Set trigonometric coefficients for forward transform
1420 if (sign == -1) {
1421 w_real = wavetable[twiddle].real();
1422 w_imag = wavetable[twiddle].imag();
1423 }
1424 // ... otherwise set trigonometric coefficients for backward
1425 // tranform: w -> conjugate(w)
1426 else {
1427 w_real = wavetable[twiddle].real();
1428 w_imag = -wavetable[twiddle].imag();
1429 }
1430 }
1431
1432 // Loop over ...
1433 for (int i = 0; i < m; ++i) {
1434
1435 //
1436 const double xp_real = (out+ostride*(i + e1 * m))->real();
1437 const double xp_imag = (out+ostride*(i + e1 * m))->imag();
1438 const double xm_real = (out+ostride*(i + (factor-e1)*m))->real();
1439 const double xm_imag = (out+ostride*(i + (factor-e1)*m))->imag();
1440
1441 //
1442 const double ap = w_real * xp_real;
1443 const double am = w_imag * xm_imag;
1444
1445 //
1446 double sum_real = ap - am;
1447 double sumc_real = ap + am;
1448
1449 //
1450 const double bp = w_real * xp_imag;
1451 const double bm = w_imag * xm_real;
1452
1453 //
1454 double sum_imag = bp + bm;
1455 double sumc_imag = bp - bm;
1456
1457 //
1458 *(in+istride*(i+em)) = std::complex<double>((in+istride*(i+em))->real() + sum_real,
1459 (in+istride*(i+em))->imag() + sum_imag);
1460 *(in+istride*(i+ecm)) = std::complex<double>((in+istride*(i+ecm))->real() + sumc_real,
1461 (in+istride*(i+ecm))->imag() + sumc_imag);
1462
1463
1464 } // endfor: i
1465
1466 // Increment
1467 idx += idx_step ;
1468 idx %= factor * q ;
1469
1470 } // endfor: e1
1471
1472 } // endfor: e
1473
1474 // k = 0
1475 for (int k1 = 0; k1 < p_1; ++k1) {
1476 *(out+ostride*k1) = *(in+istride*k1);
1477 }
1478
1479 // k > 0
1480 for (int e1 = 1; e1 < factor; ++e1) {
1481 for (int k1 = 0; k1 < p_1; ++k1) {
1482 *(out+ostride*(k1 + e1 * p_1)) = *(in+istride*(k1 + e1 * m));
1483 }
1484 }
1485
1486 // e = 0
1487 for (int k = 1, i = p_1, j = product; k < q; ++k, j += jump) {
1488 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
1489 *(out+ostride*j) = *(in+istride*i);
1490 }
1491 }
1492
1493 // e > 0
1494 for (int k = 1, i = p_1, j = product; k < q; ++k, j += jump) {
1495
1496 for (int k1 = 0; k1 < p_1; ++k1, ++i, ++j) {
1497
1498 for (int e1 = 1; e1 < factor; ++e1) {
1499
1500 // Get x
1501 std::complex<double> x = *(in+istride*(i + e1 * m));
1502
1503 // Compute index
1504 int twiddle = index + (e1-1)*q + k-1;
1505
1506 // Get w
1507 std::complex<double> w = (sign == -1)
1508 ? wavetable[twiddle]
1509 : std::conj(wavetable[twiddle]);
1510
1511 // Compute out = w * x
1512 *(out+ostride*(j + e1 * p_1)) = w * x;
1513
1514 } // endfor: e1
1515
1516 } // endfor: k1
1517
1518 } // endfor: k
1519
1520 // Return
1521 return;
1522}
1523
1524
1525/***********************************************************************//**
1526 * @brief Extract coefficients from wavetable
1527 *
1528 * @param[in] index Start index of trigonometric coefficients.
1529 * @param[in] k ...
1530 * @param[in] q ...
1531 * @param[in] n Number of coefficients.
1532 * @param[in] sign Forward (-1) or backward (+1) transformation.
1533 * @param[in] wavetable Trigonometric coefficients.
1534 ***************************************************************************/
1535std::vector<std::complex<double> > GFft::get_w(const GFftWavetable& wavetable,
1536 const int& index,
1537 const int& k,
1538 const int& q,
1539 const int& n,
1540 const int& sign) const
1541{
1542 // Allocate w vectors
1543 std::vector<std::complex<double> > w(n);
1544
1545 // Set trigonometric coefficients for k=0 since they are not stored in the
1546 // wavetable object
1547 if (k == 0) {
1548 w.assign(n, std::complex<double>(1.0, 0.0));
1549 }
1550
1551 // ... otherwise use coefficients stored in wavetable
1552 else {
1553
1554 // Compute indices
1555 int twiddle = index + k - 1;
1556
1557 // If sign==-1 then set trigonometric coefficients for forward transform
1558 if (sign == -1) {
1559 for (int i = 0; i < n; ++i, twiddle += q) {
1560 w[i] = wavetable[twiddle];
1561 }
1562 }
1563
1564 // ... otherwise set trigonometric coefficients for backward tranform
1565 // w -> conjugate(w)
1566 else {
1567 for (int i = 0; i < n; ++i, twiddle += q) {
1568 w[i] = std::conj(wavetable[twiddle]);
1569 }
1570 }
1571
1572 } // endelse: k != 0
1573
1574 // Return w vectors
1575 return w;
1576}
Exception handler interface definition.
#define G_OP_DIV
Definition GFft.cpp:41
#define G_OP_MUL
Definition GFft.cpp:40
#define G_BACKWARD
Definition GFft.cpp:43
#define G_FORWARD
Definition GFft.cpp:42
Fast Fourier transformation class interface definition.
Mathematical function definitions.
#define G_OP_ADD
Definition GMatrix.cpp:43
#define G_OP_SUB
Definition GMatrix.cpp:44
GNdarray sign(const GNdarray &array)
Computes sign of array elements.
N-dimensional array class interface definition.
Gammalib tools definition.
GChatter
Definition GTypemaps.hpp:33
@ SILENT
Definition GTypemaps.hpp:34
@ VERBOSE
Definition GTypemaps.hpp:38
double norm(const GVector &vector)
Computes vector norm.
Definition GVector.cpp:864
Lookup table class for Fast Fourier Transformation.
int factor(const int &index) const
Return factorisation factor.
const int & index(const int &factor) const
Return start index for a given factor.
int factors(void) const
Return number of factorisation factors.
Fast Fourier Transformation class.
Definition GFft.hpp:57
GFft operator-(void) const
Unary minus operator.
Definition GFft.cpp:252
GFft & operator*=(const GFft &fft)
Unary multiplication operator.
Definition GFft.cpp:206
void forward(const GNdarray &array)
Forward Fast Fourier Transform.
Definition GFft.cpp:328
int m_size
Size of data array.
Definition GFft.hpp:180
void factorn(std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &factor, const int &product, const int &n, const int &index)
Compute FFT for arbitrary factor.
Definition GFft.cpp:1343
void copy_members(const GFft &fft)
Copy class members.
Definition GFft.cpp:558
std::vector< int > m_shape
Array dimensions.
Definition GFft.hpp:182
GNdarray backward(void) const
Backward Fast Fourier Transform.
Definition GFft.cpp:401
std::complex< double > timesi(const std::complex< double > &value) const
Return complex value times i.
Definition GFft.hpp:397
void factor7(const std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &product, const int &n, const int &index)
Compute FFT for factor 7.
Definition GFft.cpp:1214
void require_same_shape(const std::string &method, const GFft &fft) const
Throw exception if FFT shapes differ.
Definition GFft.cpp:676
void set_data(const GNdarray &array)
Set data from n-dimensional array.
Definition GFft.cpp:603
virtual ~GFft(void)
Destructor.
Definition GFft.cpp:106
GFft(void)
Void constructor.
Definition GFft.cpp:55
void factor3(const std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &product, const int &n, const int &index)
Compute FFT for factor 3.
Definition GFft.cpp:885
std::vector< std::complex< double > > get_w(const GFftWavetable &wavetable, const int &index, const int &k, const int &q, const int &n, const int &sign) const
Extract coefficients from wavetable.
Definition GFft.cpp:1535
std::vector< GFftWavetable > m_wavetable
Trigonometric coefficients.
Definition GFft.hpp:184
GFft & operator+=(const GFft &fft)
Unary addition operator.
Definition GFft.cpp:158
void factor2(const std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &product, const int &n, const int &index)
Compute FFT for factor 2.
Definition GFft.cpp:825
std::string print(const GChatter &chatter=NORMAL) const
Print Fast Fourier Transform information.
Definition GFft.cpp:485
GFft & operator-=(const GFft &fft)
Unary subtraction operator.
Definition GFft.cpp:182
const std::vector< int > & shape(void) const
Return shape of array.
Definition GFft.hpp:296
int dim(void) const
Return dimension of Fast Fourier Transformation.
Definition GFft.hpp:268
void transform(std::complex< double > *data, const int &stride, const int &n, const GFftWavetable &wavetable, const bool &forward=true)
Perform Fast Fourier Transform.
Definition GFft.cpp:718
std::vector< int > m_strides
Steps in each dimension.
Definition GFft.hpp:183
void factor6(const std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &product, const int &n, const int &index)
Compute FFT for factor 6.
Definition GFft.cpp:1118
std::complex< double > * m_data
Pointer on array data.
Definition GFft.hpp:181
void clear(void)
Clear Fast Fourier Transform.
Definition GFft.cpp:276
int size(void) const
Return number of elements in array.
Definition GFft.hpp:282
void free_members(void)
Delete class members.
Definition GFft.cpp:582
bool has_same_shape(const GFft &fft) const
Check if FFT has the same shape.
Definition GFft.cpp:646
GFft & operator/=(const GFft &fft)
Unary division operator.
Definition GFft.cpp:230
void factor5(const std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &product, const int &n, const int &index)
Compute FFT for factor 5.
Definition GFft.cpp:1029
GFft & operator=(const GFft &fft)
Assignment operator.
Definition GFft.cpp:128
void factor4(const std::complex< double > *in, const int &istride, std::complex< double > *out, const int &ostride, const GFftWavetable &wavetable, const int &sign, const int &product, const int &n, const int &index)
Compute FFT for factor 4.
Definition GFft.cpp:957
void init_members(void)
Initialise class members.
Definition GFft.cpp:539
GFft * clone(void) const
Clone Fast Fourier Transform.
Definition GFft.cpp:294
N-dimensional array class.
Definition GNdarray.hpp:44
int size(void) const
Return number of elements in array.
Definition GNdarray.hpp:293
int dim(void) const
Return dimension of array.
Definition GNdarray.hpp:279
const std::vector< int > & shape(void) const
Return shape of array.
Definition GNdarray.hpp:307
const std::vector< int > & strides(void) const
Return strides of array.
Definition GNdarray.hpp:321
std::string parformat(const std::string &s, const int &indent=0)
Convert string in parameter format.
Definition GTools.cpp:1143
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
Definition GTools.cpp:489
const double twopi
Definition GMath.hpp:36