NAG CL Interface
C06 (Sum)
Fourier Transforms

Settings help

CL Name Style:


1 Scope of the Chapter

This chapter is concerned with the following tasks:
  1. (a)calculating the discrete Fourier transform of a sequence of real or complex data values;
  2. (b)calculating the discrete convolution or the discrete correlation of two sequences of real or complex data values using discrete Fourier transforms;
  3. (c)calculating the fast Gauss transform approximation to the discrete Gauss transform;
  4. (d)direct summation of orthogonal series;

2 Background to the Problems

2.1 Discrete Fourier Transforms

2.1.1 Complex transforms

Most of the functions in this chapter calculate the finite discrete Fourier transform (DFT) of a sequence of n complex numbers zj , for j=0,1,,n-1. The direct transform is defined by
z^k = 1n j=0 n-1 zj exp(-i 2πjk n ) (1)
for k=0,1,,n-1 . Note that equation (1) makes sense for all integral k and with this extension z^k is periodic with period n, i.e., z^k = z^ k±n , and in particular z^ -k = z^ n-k . Note also that the scale-factor of 1n may be omitted in the definition of the DFT, and replaced by 1n in the definition of the inverse.
If we write zj = xj + i yj and z^k = ak + i bk , then the definition of z^k may be written in terms of sines and cosines as
ak = 1n j=0 n-1 (xjcos( 2πjk n )+yjsin( 2πjk n ))  
bk = 1n j= 0 n- 1 (yjcos( 2πjk n )-xjsin( 2πjk n )) .  
The original data values zj may conversely be recovered from the transform z^k by an inverse discrete Fourier transform:
zj = 1n k=0 n-1 z^k exp(+i 2πjk n ) (2)
for j=0,1,,n-1 . If we take the complex conjugate of (2), we find that the sequence z¯j is the DFT of the sequence z^¯ k . Hence the inverse DFT of the sequence z^k may be obtained by taking the complex conjugates of the z^k ; performing a DFT, and taking the complex conjugates of the result. (Note that the terms forward transform and backward transform are also used to mean the direct and inverse transforms respectively.)
The definition (1) of a one-dimensional transform can easily be extended to multidimensional transforms. For example, in two dimensions we have
z^ k1 k2 = 1 n1 n2 j1=0 n1-1 j2=0 n2-1 z j1 j2 exp(-i 2 π j1 k1 n1 ) exp(-i 2 π j2 k2 n2 ) . (3)
Note:  definitions of the discrete Fourier transform vary. Sometimes (2) is used as the definition of the DFT, and (1) as the definition of the inverse.

2.1.2 Real transforms

If the original sequence is purely real valued, i.e., zj=xj , then
z^k = ak + i bk = 1n j=0 n-1 xj exp(-i 2πjk n )  
and z^ n-k is the complex conjugate of z^k . Thus the DFT of a real sequence is a particular type of complex sequence, called a Hermitian sequence, or half-complex or conjugate symmetric, with the properties
a n-k = ak b n-k = -bk b0 = 0  
and, if n is even, b n/2 = 0 .
Thus a Hermitian sequence of n complex data values can be represented by only n, rather than 2n , independent real values. This can obviously lead to economies in storage, with two schemes being used in this chapter. In the first (deprecated) scheme, which will be referred to as the real storage format for Hermitian sequences, the real parts ak for 0 k n/2 are stored in normal order in the first n/2 + 1 locations of an array x of length n; the corresponding nonzero imaginary parts are stored in reverse order in the remaining locations of x. To clarify, the following two tables illustrate the storage of the real and imaginary parts of z^k for the two cases: n even and n odd.
If n is even then the sequence has two purely real elements and is stored as follows:
Index of x 0 1 2 n/2 n-2 n-1
Sequence a0 a1 + ib1 a2 + ib2 a n/2 a2 - ib2 a1 - ib1
Stored values a0 a1 a2 a n/2 b2 b1
x[k] = ak , for ​ k= 0, 1, , n/2 , and x[n-k] = bk , for ​ k= 1, 2, , n/2-1 .  
If n is odd then the sequence has one purely real element and, letting n=2s+1 , is stored as follows:
Index of x 0 1 2 s s+1 n-2 n-1
Sequence a0 a1 + ib1 a2 + ib2 as + ibs as - ibs a2 - ib2 a1 - ib1
Stored values a0 a1 a2 as bs b2 b1
x[k] = ak , for ​ k= 0, 1, , s , and x[n-k] = bk , for ​ k= 1, 2, , s .  
The second (recommended) storage scheme, referred to in this chapter as the complex storage format for Hermitian sequences, stores the real and imaginary parts ak , bk , for 0 k n/2 , in consecutive locations of an array x of length n+2 . The following two tables illustrate the storage of the real and imaginary parts of z^k for the two cases: n even and n odd.
If n is even then the sequence has two purely real elements and is stored as follows:
Index of x 0 1 2 3 n-2 n-1 n n+1
Stored values a0 b0=0 a1 b1 a n/2-1 b n/2-1 a n/2 b n/2 = 0
x[2×k-1] = ak , for ​ k= 0, 1, , n/2 , and x[2×k+1-1] = bk , for ​ k= 0, 1, , n/2 .  
If n is odd then the sequence has one purely real element and, letting n=2s+1 , is stored as follows:
Index of x 0 1 2 3 n-2 n-1 n n+1
Stored values a0 b0=0 a1 b1 b s-1 as bs 0
x[2×k-1] = ak , for ​ k= 0, 1, , s , and x[2×k+1-1] = bk , for ​ k= 0, 1, , s .  
Also, given a Hermitian sequence, the inverse (or backward) discrete transform produces a real sequence. That is,
xj = 1n (a0+2 k=1 n/2-1 (akcos( 2πjk n )-bksin( 2πjk n ))+ a n/2 )  
where a n/2 = 0 if n is odd.
For real data that is two-dimensional or higher, the symmetry in the transform persists for the leading dimension only. So, using the notation of equation (3) for the complex two-dimensional discrete transform, we have that z^k1k2 is the complex conjugate of z^(n1-k1)(n2-k2). It is more convenient for transformed data of two or more dimensions to be stored as a complex sequence of length (n1/2+1)×n2××nd where d is the number of dimensions. The inverse discrete Fourier transform operating on such a complex sequence (Hermitian in the leading dimension) returns a real array of full dimension (n1×n2××nd).

2.1.3 Real symmetric transforms

In many applications the sequence xj will not only be real, but may also possess additional symmetries which we may exploit to reduce further the computing time and storage requirements. For example, if the sequence xj is odd, (xj= -x n-j ) , then the discrete Fourier transform of xj contains only sine terms. Rather than compute the transform of an odd sequence, we define the sine transform of a real sequence by
x^k = 2n j=1 n-1 xj sin(πjkn) ,  
which could have been computed using the Fourier transform of a real odd sequence of length 2n . In this case the xj are arbitrary, and the symmetry only becomes apparent when the sequence is extended. Similarly we define the cosine transform of a real sequence by
x^k = 2n (12x0+ j=1 n-1 xjcos(πjkn)+12(-1)kxn)  
which could have been computed using the Fourier transform of a real even sequence of length 2n .
In addition to these ‘half-wave’ symmetries described above, sequences arise in practice with ‘quarter-wave’ symmetries. We define the quarter-wave sine transform by
x^k = 1n ( j=1 n-1 xjsin( π j (2k-1) 2n )+12(-1) k-1xn)  
which could have been computed using the Fourier transform of a real sequence of length 4n of the form
(0,x1,,xn,xn-1 ,,x1,0,-x1,,-xn, -x n-1 ,, -x 1 ) .  
Similarly we may define the quarter-wave cosine transform by
x^k = 1n (12x0+j= 1 n- 1xjcos( π j (2k-1) 2n ))  
which could have been computed using the Fourier transform of a real sequence of length 4n of the form
(x0,x1,, x n-1 ,0, -x n-1 ,,-x0,-x1,, -x n-1 ,0, x n-1 ,,x1) .  

2.1.4 Fourier integral transforms

The usual application of the discrete Fourier transform is that of obtaining an approximation of the Fourier integral transform
F (s) = - f (t) exp(-i2πst) dt  
when f(t) is negligible outside some region (0,c) . Dividing the region into n equal intervals we have
F (s) cn j=0 n-1 fj exp( -i 2 π s j c n )  
and so
Fk cn j= 0 n- 1 fj exp( -i 2 π jk n )  
for k=0,1,,n-1 , where fj = f (jc/n) and Fk = F (k/c) .
Hence the discrete Fourier transform gives an approximation to the Fourier integral transform in the region s=0 to s=n/c .
If the function f(t) is defined over some more general interval (a,b) , then the integral transform can still be approximated by the discrete transform provided a shift is applied to move the point a to the origin.

2.1.5 Convolutions and correlations

One of the most important applications of the discrete Fourier transform is to the computation of the discrete convolution or correlation of two vectors x and y defined (as in Brigham (1974)) by (Here x and y are assumed to be periodic with period n.)
Under certain circumstances (see Brigham (1974)) these can be used as approximations to the convolution or correlation integrals defined by
z (s) = - x (t) y (s-t) dt  
and
w (s) = - x¯ (t) y (s+t) dt ,   - < s < .  
For more general advice on the use of Fourier transforms, see Hamming (1962); more detailed information on the fast Fourier transform algorithm can be found in Gentleman and Sande (1966) and Brigham (1974).

2.1.6 Applications to solving partial differential equations (PDEs)

A further application of the fast Fourier transform, and in particular of the Fourier transforms of symmetric sequences, is in the solution of elliptic PDEs. If an equation is discretized using finite differences, then it is possible to reduce the problem of solving the resulting large system of linear equations to that of solving a number of tridiagonal systems of linear equations. This is accomplished by uncoupling the equations using Fourier transforms, where the nature of the boundary conditions determines the choice of transforms – see Section 3.3. Full details of the Fourier method for the solution of PDEs may be found in Swarztrauber (1977) and Swarztrauber (1984).

2.2 Fast Gauss Transform

Gauss transforms have applications in areas including statistics, machine learning, and numerical solution of the heat equation. The discrete Gauss transform (DGT), G(y), evaluated at a set of target points y(j), for j=1,2,,md, is defined as:
G (yj) = i=1 n qi e - yj-xi 2 2 / h i 2 ,  j=1,,m  
where xi, for i=1,2,,nd, are the Gaussian source points, qi, for i=1,2,,n+, are the source weights and hi, for i=1,2,,n+, are the source standard deviations (alternatively source scales or source bandwidths).
The fast Gauss transform (FGT) algorithm presented in Raykar and Duraiswami (2005) approximates the DGT by using two Taylor series and clustering of the source points.

2.3 Direct Summation of Orthogonal Series

For any series of functions ϕi which satisfy a recurrence
ϕr+1 (x) + αr (x) ϕr (x) + βr (x) ϕr-1 (x) =0  
the sum
r= 0 n ar ϕr (x)  
is given by
r=0 n ar ϕr (x) = b0 (x) ϕ0 (x) + b1 (x) (ϕ1(x)+α0(x)ϕ0(x))  
where
br (x) + αr (x) br+ 1 (x) + βr+ 1 (x) br+ 2 (x) = ar b n+ 1 (x) = b n+ 2 (x) = 0 .  
This may be used to compute the sum of the series. For further reading, see Hamming (1962).

3 Recommendations on Choice and Use of Available Functions

The fast Fourier transform algorithm ceases to be ‘fast’ if applied to values of n which cannot be expressed as a product of small prime factors. All the FFT functions in this chapter are particularly efficient if the only prime factors of n are 2, 3 or 5.

3.1 One-dimensional Fourier Transforms

The choice of function is determined first of all by whether the data values constitute a real, Hermitian or general complex sequence. It is wasteful of time and storage to use an inappropriate function.

3.1.1 Real and Hermitian data

c06pac transforms a single sequence of real data onto (and in-place) a representation of the transformed Hermitian sequence using the complex storage scheme described in Section 2.1.2. c06pac also performs the inverse transform using the representation of Hermitian data and transforming back to a real data sequence.
Alternatively, the two-dimensional function c06pvc can be used (on setting the second dimension to 1) to transform a sequence of real data onto an Hermitian sequence whose first half is stored in a separate Complex array. The second half need not be stored since these are the complex conjugate of the first half in reverse order. c06pwc performs the inverse operation, transforming the the Hermitian sequence (half-)stored in a Complex array onto a separate real array.
Multiple sequences of real data can also be transformed to and from a transformed Hermitian sequence (using the complex storage scheme). c06ppc assumes that the original real data is stored such that all the corresponding (j-th say) elements of the different streams are stored contiguously. c06pqc assumes that the elements of each sequence are stored contiguously; this is the recommended form of storage since it normally results in more efficient computation.

3.1.2 Complex data

c06pcc transforms a single complex sequence in-place; it also performs the inverse transform. c06psc transforms multiple complex sequences, each stored sequentially; it also performs the inverse transform on multiple complex sequences. This function is designed to perform several transforms in a single call, all with the same value of n.
If extensive use is to be made of these functions and you are concerned about efficiency, you are advised to conduct your own timing tests.

3.2 Half- and Quarter-wave Transforms

Four functions are provided for computing fast Fourier transforms (FFTs) of real symmetric sequences. c06rec computes multiple Fourier sine transforms, c06rfc computes multiple Fourier cosine transforms, c06rgc computes multiple quarter-wave Fourier sine transforms, and c06rhc computes multiple quarter-wave Fourier cosine transforms.

3.3 Application to Elliptic Partial Differential Equations

As described in Section 2.1.6, Fourier transforms may be used in the solution of elliptic PDEs.
c06rec may be used to solve equations where the solution is specified along the boundary.
c06rfc may be used to solve equations where the derivative of the solution is specified along the boundary.
c06rgc may be used to solve equations where the solution is specified on the lower boundary, and the derivative of the solution is specified on the upper boundary.
c06rhc may be used to solve equations where the derivative of the solution is specified on the lower boundary, and the solution is specified on the upper boundary.
For equations with periodic boundary conditions the full-range Fourier transforms computed by c06pac are appropriate.

3.4 Multidimensional Fourier Transforms

The following functions compute multidimensional discrete Fourier transforms of real, Hermitian and complex data stored in Complex arrays:
double Hermitian Complex
2 dimensions c06pvc c06pwc c06puc
3 dimensions c06pyc c06pzc c06pxc
any number of dimensions c06pjc
The Hermitian data, either transformed from or being transformed to real data, is compacted (due to symmetry) along its first dimension when stored in Complex arrays; thus approximately half the full Hermitian data is stored.
c06puc and c06pxc should be used in preference to c06pjc for two- and three-dimensional transforms, as they are easier to use and are likely to be more efficient.
The transform of multidimensional real data is stored as a complex sequence that is Hermitian in its leading dimension. The inverse transform takes such a complex sequence and computes the real transformed sequence. Consequently, separate functions are provided for performing forward and inverse transforms.
c06pvc performs the forward two-dimensionsal transform while c06pwc performs the inverse of this transform.
c06pyc performs the forward three-dimensional transform while c06pzc performs the inverse of this transform.
The complex sequences computed by c06pvc and c06pyc contain roughly half of the Fourier coefficients; the remainder can be reconstructed by conjugation of those computed. For example, the Fourier coefficients of the two-dimensional transform z^(n1-k1)k2 are the complex conjugate of z^k1k2 for k1=0,1,,n1/2, and k2=0,1,,n2-1.

3.5 Convolution and Correlation

c06fkc computes either the discrete convolution or the discrete correlation of two real vectors.

3.6 Fast Gauss Transform

The only function available is c06sac. If the dimensionality of the data is low or the number of source and target points is small, however, it may be more efficient to evaluate the discrete Gauss transform directly.

3.7 Direct Summation of Orthogonal Series

The only function available is c06dcc, which sums a finite Chebyshev series
j=0 n cj Tj (x) ,   j=0 n cj T 2j (x)   or   j=0 n cj T 2j+1 (x)  
depending on the choice of argument.

4 Decision Trees

Tree 1: Fourier Transform of Discrete Complex Data

Is the data one-dimensional?   Multiple vectors?   c06psc
yesyes
  no   no
c06pcc
Is the data two-dimensional?   c06puc
yes
  no
Is the data three-dimensional?   c06pxc
yes
  no
Transform on one dimension only?   c06pfc
yes
  no
Transform on all dimensions?   c06pjc
yes

Tree 2: Fourier Transform of Real Data or Data in Complex Hermitian Form Resulting from the Transform of Real Data

Quarter-wave sine (inverse) transform?   c06rgc
yes
  no
Quarter-wave cosine (inverse) transform?   c06rhc
yes
  no
Sine (inverse) transform?   c06rec
yes
  no
Cosine (inverse) transform?   c06rfc
yes
  no
Is the data three-dimensional?   Forward transform on real data?   c06pyc
yesyes
  no   no
Inverse transform on Hermitian data?   c06pzc
yes
Is the data two-dimensional?   Forward transform on real data?   c06pvc
yesyes
  no   no
Inverse transform on Hermitian data?   c06pwc
yes
Is the data multi one-dimensional?   Sequence elements stored by stride (number of sequences)?   c06ppc
yesyes
  no   no
Sequence elements stored contiguously?   c06pqc
yes
c06pac

5 Functionality Index

Convolution or Correlation,  
real vectors,  
time-saving   c06fkc
Discrete Fourier Transform,  
multidimensional,  
complex sequence,  
complex storage   c06pjc
multiple half- and quarter-wave transforms,  
Fourier cosine transforms, simple use   c06rfc
Fourier sine transforms, simple use   c06rec
quarter-wave cosine transforms, simple use   c06rhc
quarter-wave sine transforms, simple use   c06rgc
one-dimensional,  
multiple transforms,  
complex sequence,  
complex storage, contiguous sequence elements   c06psc
Hermitian/real sequence,  
complex storage, contiguous sequence elements   c06pqc
complex storage, strided sequence elements   c06ppc
multi-variable,  
complex sequence,  
complex storage   c06pfc
single transforms,  
complex sequence,  
time-saving,  
complex storage   c06pcc
Hermitian/real sequence,  
time-saving,  
complex storage   c06pac
three-dimensional,  
complex sequence,  
complex storage   c06pxc
Hermitian/real sequence,  
complex-to-real   c06pzc
real-to-complex   c06pyc
two-dimensional,  
complex sequence,  
complex storage   c06puc
Hermitian/real sequence,  
complex-to-real   c06pwc
real-to-complex   c06pvc
Fast Gauss Transform   c06sac
Summation of Chebyshev series   c06dcc

6 Auxiliary Functions Associated with Library Function Arguments

None.

7 Withdrawn or Deprecated Functions

The following lists all those functions that have been withdrawn since Mark 23 of the Library or are in the Library, but deprecated.
Function Status Replacement Function(s)
c06eac Withdrawn at Mark 26 c06pac
c06ebc Withdrawn at Mark 26 c06pac
c06ecc Withdrawn at Mark 26 c06pcc
c06ekc Withdrawn at Mark 26 c06fkc
c06fpc To be withdrawn at Mark 29 c06pqc
c06fqc To be withdrawn at Mark 29 c06pqc
c06frc Withdrawn at Mark 26 c06psc
c06fuc Withdrawn at Mark 26 c06puc
c06gbc Withdrawn at Mark 26 No replacement required
c06gcc Withdrawn at Mark 26 No replacement required
c06gqc To be withdrawn at Mark 29 No replacement required
c06gsc To be withdrawn at Mark 29 No replacement required
c06gzc To be withdrawn at Mark 29 No replacement required
c06hac Withdrawn at Mark 26 c06rec
c06hbc Withdrawn at Mark 26 c06rfc
c06hcc Withdrawn at Mark 26 c06rgc
c06hdc Withdrawn at Mark 26 c06rhc

8 References

Brigham E O (1974) The Fast Fourier Transform Prentice–Hall
Davies S B and Martin B (1979) Numerical inversion of the Laplace transform: A survey and comparison of methods J. Comput. Phys. 33 1–32
Fox L and Parker I B (1968) Chebyshev Polynomials in Numerical Analysis Oxford University Press
Gentleman W S and Sande G (1966) Fast Fourier transforms for fun and profit Proc. Joint Computer Conference, AFIPS 29 563–578
Hamming R W (1962) Numerical Methods for Scientists and Engineers McGraw–Hill
Raykar V C and Duraiswami R (2005) Improved Fast Gauss Transform With Variable Source Scales University of Maryland Technical Report CS-TR-4727/UMIACS-TR-2005-34
Shanks D (1955) Nonlinear transformations of divergent and slowly convergent sequences J. Math. Phys. 34 1–42
Swarztrauber P N (1977) The methods of cyclic reduction, Fourier analysis and the FACR algorithm for the discrete solution of Poisson's equation on a rectangle SIAM Rev. 19(3) 490–501
Swarztrauber P N (1984) Fast Poisson solvers Studies in Numerical Analysis (ed G H Golub) Mathematical Association of America
Swarztrauber P N (1986) Symmetric FFT's Math. Comput. 47(175) 323–346
Wynn P (1956) On a device for computing the em(Sn) transformation Math. Tables Aids Comput. 10 91–96