casa  $Rev:20696$
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
SCSL.h
Go to the documentation of this file.
00001 //# extern_fft.h: C++ wrapper functions for FORTRAN FFT code
00002 //# Copyright (C) 1993,1994,1995,1997,1999,2000
00003 //# Associated Universities, Inc. Washington DC, USA.
00004 //#
00005 //# This library is free software; you can redistribute it and/or modify it
00006 //# under the terms of the GNU Library General Public License as published by
00007 //# the Free Software Foundation; either version 2 of the License, or (at your
00008 //# option) any later version.
00009 //#
00010 //# This library is distributed in the hope that it will be useful, but WITHOUT
00011 //# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00012 //# FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Library General Public
00013 //# License for more details.
00014 //#
00015 //# You should have received a copy of the GNU Library General Public License
00016 //# along with this library; if not, write to the Free Software Foundation,
00017 //# Inc., 675 Massachusetts Ave, Cambridge, MA 02139, USA.
00018 //#
00019 //# Correspondence concerning AIPS++ should be addressed as follows:
00020 //#        Internet email: aips2-request@nrao.edu.
00021 //#        Postal address: AIPS++ Project Office
00022 //#                        National Radio Astronomy Observatory
00023 //#                        520 Edgemont Road
00024 //#                        Charlottesville, VA 22903-2475 USA
00025 //#
00026 //# $Id: SCSL.h 18093 2004-11-30 17:51:10Z ddebonis $
00027 
00028 #ifndef SCIMATH_SCSL_H
00029 #define SCIMATH_SCSL_H
00030 
00031 #include <casa/aips.h>
00032 #include <casa/BasicSL/Complex.h>
00033 
00034 
00035 namespace casa { //# NAMESPACE CASA - BEGIN
00036 
00037 // <summary>C++ Interface to the Sgi/Cray Scientific Library (SCSL)</summary>
00038 // <synopsis>
00039 // These are C++ wrapper functions for the transform routines in the SGI/Cray 
00040 // Scientific Library (SCSL). The purpose of these definitions is to overload
00041 // the functions so that C++ users can access the functions in SCSL with
00042 // identical function names.
00043 //
00044 // <note role=warning> 
00045 // Currently, the SCSL is available only on SGI machines.
00046 // </note>
00047 // </synopsis>
00048 
00049 class SCSL
00050 {
00051 public:
00052 // These routines compute the Fast Fourier Transform (FFT) of the complex
00053 // vector x, and store the result in vector y.  <src>ccfft</src> does the
00054 // complex-to-complex transform and <src>zzfft</src> does the same for double
00055 // precision arrays.
00056 //
00057 // In FFT applications, it is customary to use zero-based subscripts; the
00058 // formulas are simpler that way.  Suppose that the arrays are
00059 // dimensioned as follows:
00060 //
00061 // <srcblock>
00062 //      COMPLEX X(0:N-1), Y(0:N-1)
00063 // </srcblock>
00064 //
00065 // The output array is the FFT of the input array, using the following
00066 // formula for the FFT:
00067 //
00068 // <srcblock>
00069 //                     n-1
00070 //      Y(k) = scale * Sum [ X(j)*w**(isign*j*k) ]    for k = 0, ..., n-1
00071 //                     j=0
00072 //
00073 //      where:
00074 //      w = exp(2*pi*i/n),
00075 //      i = + sqrt(-1),
00076 //      pi = 3.14159...,
00077 //      isign = +1 or -1
00078 // </srcblock>
00079 //
00080 // Different authors use different conventions for which of the
00081 // transforms, isign = +1 or isign = -1, is the forward or inverse
00082 // transform, and what the scale factor should be in either case.  You
00083 // can make this routine compute any of the various possible definitions,
00084 // however, by choosing the appropriate values for isign and scale.
00085 //
00086 // The relevant fact from FFT theory is this:  If you take the FFT with
00087 // any particular values of isign and scale, the mathematical inverse
00088 // function is computed by taking the FFT with -isign and 1/(n*scale).
00089 // In particular, if you use isign = +1 and scale = 1.0 you can compute
00090 // the inverse FFT by using isign = -1 and scale = 1.0/n.
00091 //
00092 // The output array may be the same as the input array.
00093 //
00094 // <h3>Initialization</h3>
00095 // The table array stores the trigonometric tables used in calculation of
00096 // the FFT.  You must initialize table by calling the routine with isign
00097 // = 0 prior to doing the transforms.  If the value of the problem size,
00098 // n, does not change, table does not have to be reinitialized.
00099 //
00100 // <h3>Dimensions</h3>
00101 //  In the preceding description, it is assumed that array subscripts were
00102 // zero-based, as is customary in FFT applications.  Thus, the input and
00103 // output arrays are declared as follows:
00104 //
00105 // <srcblock>
00106 //      COMPLEX X(0:N-1)
00107 //      COMPLEX Y(0:N-1)
00108 // </srcblock>
00109 //
00110 // However, if you prefer to use the more customary FORTRAN style with
00111 // subscripts starting at 1 you do not have to change the calling
00112 // sequence, as in the following (assuming N > 0):
00113 //
00114 // <srcblock>
00115 //      COMPLEX X(N)
00116 //      COMPLEX Y(N)
00117 // </srcblock>
00118 //
00119 // <example>
00120 // These examples use the table and workspace sizes appropriate to the
00121 // Origin series.
00122 //
00123 // Example 1:  Initialize the complex array table in preparation for
00124 // doing an FFT of size 1024.  Only the isign, n, and table arguments are
00125 // used in this case.  You can use dummy arguments or zeros for the other
00126 // arguments in the subroutine call.
00127 //
00128 // <srcblock> 
00129 //      REAL TABLE(30 + 2048)
00130 //      CALL CCFFT(0, 1024, 0.0, DUMMY, DUMMY, TABLE, DUMMY, 0)
00131 // </srcblock> 
00132 //
00133 // Example 2:  x and y are complex arrays of dimension (0:1023).  Take
00134 // the FFT of x and store the results in y.  Before taking the FFT,
00135 // initialize the table array, as in example 1.
00136 //
00137 // <srcblock> 
00138 //      COMPLEX X(0:1023), Y(0:1023)
00139 //      REAL TABLE(30 + 2048)
00140 //      REAL WORK(2048)
00141 //      CALL CCFFT(0, 1024, 1.0, X, Y, TABLE, WORK, 0)
00142 //      CALL CCFFT(1, 1024, 1.0, X, Y, TABLE, WORK, 0)
00143 // </srcblock> 
00144 //
00145 // Example 3:  Using the same x and y as in example 2, take the inverse
00146 // FFT of y and store it back in x.  The scale factor 1/1024 is used.
00147 // Assume that the table array is already initialized.
00148 //
00149 // <srcblock> 
00150 //      CALL CCFFT(-1, 1024, 1.0/1024.0, Y, X, TABLE, WORK, 0)
00151 // </srcblock> 
00152 //
00153 // Example 4:  Perform the same computation as in example 2, but assume
00154 // that the lower bound of each array is 1, rather than 0.  No change was
00155 // needed in the subroutine calls.
00156 //
00157 // <srcblock> 
00158 //      COMPLEX X(1024), Y(1024)
00159 //      CALL CCFFT(0, 1024, 1.0, X, Y, TABLE, WORK, 0)
00160 //      CALL CCFFT(1, 1024, 1.0, X, Y, TABLE, WORK, 0)
00161 // </srcblock> 
00162 //
00163 // Example 5:  Do the same computation as in example 4, but put the
00164 // output back in array x to save storage space.  Assume that table is
00165 // already initialized.
00166 //
00167 // <srcblock> 
00168 //      COMPLEX X(1024)
00169 //      CALL CCFFT(1, 1024, 1.0, X, X, TABLE, WORK, 0)
00170 // </srcblock> 
00171 // </example>
00172 //
00173 // Input parameters:
00174 // <dl compact>
00175 // <dt><b>isign</b>
00176 // <dd>    Integer.
00177 //         Specifies whether to initialize the table array or to do the
00178 //         forward or inverse Fourier transform, as follows:
00179 //
00180 //         If isign = 0, the routine initializes the table array and
00181 //         returns.  In this case, the only arguments used or checked
00182 //         are isign, n, and table.
00183 //
00184 //         If isign = +1 or -1, the value of isign is the sign of the
00185 //         exponent used in the FFT formula.
00186 // <dt><b>n</b>
00187 // <dd>    Integer.  Size of the transform (the number of values in
00188 //         the input array).  n >= 1.
00189 // <dt><b>scale</b>
00190 // <dd>    Scale factor.  
00191 //         <src>ccfft</src>: real.
00192 //         <src>zzfft</src>: double precision.
00193 //         Each element of the output array is multiplied by scale
00194 //         after taking the Fourier transform, as defined by the previous
00195 //         formula.
00196 // <dt><b>x</b>
00197 // <dd>    Array of dimension (0:n-1).
00198 //         <src>ccfft</src>: complex array.
00199 //         <src>zzfft</src>: double complex array. 
00200 // 
00201 //         Input array of values to be transformed.
00202 // <dt><b>isys</b>
00203 // <dd>    Integer.
00204 //         Algorithm used; value dependent on hardware system.  Currently, no
00205 //         special options are supported; therefore, you must always specify
00206 //         an isys argument as constant 0.
00207 // </dl>
00208 // Output parameters:
00209 // <dl compact>
00210 // <dt><b>y</b>
00211 // <dd>    Array of dimension (0:n-1).
00212 //         <src>ccfft</src>: complex array.
00213 //         <src>zzfft</src>: double complex array.
00214 //         Output array of transformed values.  The output array may be
00215 //         the same as the input array.  In that case, the transform is
00216 //         done in place and the input array is overwritten with the
00217 //         transformed values.
00218 // <dt><b>table</b>
00219 // <dd>    Real array; dimension 2*n+30.
00220 //
00221 //         Table of factors and trigonometric functions.
00222 //
00223 //         If isign = 0, the routine initializes table (table is output
00224 //         only).
00225 //
00226 //         If isign = +1 or -1, the values in table are assumed to be
00227 //         initialized already by a prior call with isign = 0 (table is
00228 //         input only).
00229 // <dt><b>work</b>
00230 // <dd>    Real array; dimension 2*n.
00231 //
00232 //         Work array.  This is a scratch array used for intermediate
00233 //         calculations.  Its address space must be different address
00234 //         space from that of the input and output arrays.
00235 // </dl>
00236 // <group>
00237 static void ccfft(Int isign, Int n, Float scale, Complex* x,
00238                   Complex* y, Float* table, Float* work, Int isys);
00239 static void ccfft(Int isign, Int n, Double scale, DComplex* x,
00240                   DComplex* y, Double* table, Double* work, Int isys);
00241 static void zzfft(Int isign, Int n, Double scale, DComplex* x,
00242                   DComplex* y, Double* table, Double* work, Int isys);
00243 // </group>
00244 
00245 // <src>scfft/dzfft</src> computes the FFT of the real array x, and it stores
00246 // the results in the complex array y.  <src>csfft/zdfft</src> computes the
00247 // corresponding inverse complex-to-real transform.
00248 //
00249 // It is customary in FFT applications to use zero-based subscripts; the
00250 // formulas are simpler that way.  For these routines, suppose that the
00251 // arrays are dimensioned as follows:
00252 //
00253 // <srcblock>
00254 //      REAL    X(0:n-1)
00255 //      COMPLEX Y(0:n/2)
00256 // </srcblock>
00257 //
00258 // Then the output array is the FFT of the input array, using the
00259 // following formula for the FFT:
00260 //
00261 // <srcblock> 
00262 //                     n-1
00263 //      Y(k) = scale * Sum [ X(j)*w**(isign*j*k) ]    for k = 0, ..., n/2
00264 //                     j=0
00265 //
00266 //      where:
00267 //      w = exp(2*pi*i/n),
00268 //      i = + sqrt(-1),
00269 //      pi = 3.14159...,
00270 //      isign = +1 or -1.
00271 // </srcblock> 
00272 //
00273 // Different authors use different conventions for which of the
00274 // transforms, isign = +1 or isign = -1, is the forward or inverse
00275 // transform, and what the scale factor should be in either case.  You
00276 // can make these routines compute any of the various possible
00277 // definitions, however, by choosing the appropriate values for isign and
00278 // scale.
00279 //
00280 // The relevant fact from FFT theory is this:  If you call <src>scfft</src> 
00281 // with any particular values of isign and scale, the mathematical inverse
00282 // function is computed by calling <src>csfft</src> with -isign and 
00283 // 1/(n*scale).  In particular, if you use isign = +1 and scale = 1.0 in 
00284 // <src>scfft</src> for the forward FFT, you can compute the inverse FFT by
00285 // using <src>ccfft</src> with isign = -1 and scale = 1.0/n.
00286 //
00287 // <h3>Real-to-complex FFTs</h3>
00288 // Notice in the preceding formula that there are n real input values,
00289 // and n/2 + 1 complex output values.  This property is characteristic of
00290 // real-to-complex FFTs.
00291 //
00292 // The mathematical definition of the Fourier transform takes a sequence
00293 // of n complex values and transforms it to another sequence of n complex
00294 // values.  A complex-to-complex FFT routine, such as <src>ccfft</src>, will
00295 // take n complex input values, and produce n complex output values.  In
00296 // fact, one easy way to compute a real-to-complex FFT is to store the
00297 // input data in a complex array, then call routine <src>ccfft</src> to
00298 // compute the FFT.  You get the same answer when using the <src>scfft</src>
00299 // routine.
00300 //
00301 // The reason for having a separate real-to-complex FFT routine is
00302 // efficiency.  Because the input data is real, you can make use of this
00303 // fact to save almost half of the computational work.  The theory of
00304 // Fourier transforms tells us that for real input data, you have to
00305 // compute only the first n/2 + 1 complex output values, because the
00306 // remaining values can be computed from the first half of the values by
00307 // the simple formula:
00308 //
00309 // <srcblock> 
00310 //      Y(k) = conjg(Y(n-k)) for n/2 <= k <= n-1
00311 // </srcblock>
00312 //
00313 // where the notation conjgY represents the complex conjugate of y.
00314 //
00315 // In fact, in many applications, the second half of the complex output
00316 // data is never explicitly computed or stored.  Likewise, as explained
00317 // later, only the first half of the complex data has to be supplied for
00318 // the complex-to-real FFT.
00319 //
00320 // Another implication of FFT theory is that, for real input data, the
00321 // first output value, Y(0), will always be a real number; therefore, the
00322 // imaginary part will always be 0.  If n is an even number, Y(n/2) will
00323 // also be real and thus, have zero imaginary parts.
00324 //
00325 // <h3>Complex-to-real FFTs</h3>
00326 // Consider the complex-to-real case.  The effect of the computation is
00327 // given by the preceding formula, but with X complex and Y real.
00328 //
00329 // Generally, the FFT transforms a complex sequence into a complex
00330 // sequence.  However, in a certain application we may know the output
00331 // sequence is real.  Often, this is the case because the complex input
00332 // sequence was the transform of a real sequence.  In this case, you can
00333 // save about half of the computational work.
00334 //
00335 // According to the theory of Fourier transforms, for the output
00336 // sequence, Y, to be a real sequence, the following identity on the
00337 // input sequence, X, must be true:
00338 //
00339 // <srcblock> 
00340 //      X(k) = conjg(X(n-k)) for n/2 <= k <= n-1
00341 // </srcblock> 
00342 //
00343 // And, in fact, the input values X(k) for k > n/2 need not be supplied;
00344 // they can be inferred from the first half of the input.
00345 //
00346 // Thus, in the complex-to-real routine, <src>csfft</src>, the arrays can be
00347 // dimensioned as follows:
00348 //
00349 // <srcblock> 
00350 //      COMPLEX X(0:n/2)
00351 //      REAL    Y(0:n-1)
00352 // </srcblock> 
00353 //
00354 // There are n/2 + 1 complex input values and n real output values.  Even
00355 // though only n/2 + 1 input values are supplied, the size of the
00356 // transform is still n in this case, because implicitly you are using
00357 // the FFT formula for a sequence of length n.
00358 //
00359 // Another implication of the theory is that X(0) must be a real number
00360 // (that is, it must have zero imaginary part).  Also, if n is even,
00361 // X(n/2) must also be real.  Routine <src>CSFFT</src> assumes that these
00362 // values are real; if you specify a nonzero imaginary part, it is ignored.
00363 //
00364 // <h3>Table Initialization</h3>
00365 // The table array stores the trigonometric tables used in calculation of
00366 // the FFT.  This table must be initialized by calling the routine with
00367 // isign = 0 prior to doing the transforms.  The table does not have to
00368 // be reinitialized if the value of the problem size, n, does not change.
00369 // Because <src>scfft</src> and <src>csfft</src> use the same format for 
00370 // table, either can be used to initialize it (note that CCFFT uses a
00371 // different table format).
00372 //
00373 // <h3>Dimensions</h3>
00374 // In the preceding description, it is assumed that array subscripts were
00375 // zero-based, as is customary in FFT applications.  Thus, the input and
00376 // output arrays are declared (assuming n > 0):
00377 //
00378 // <srcblock> 
00379 //      REAL    X(0:n-1)
00380 //      COMPLEX Y(0:n/2)
00381 // </srcblock> 
00382 //
00383 // No change is needed in the calling sequence; however, if you prefer
00384 // you can use the more customary Fortran style with subscripts starting
00385 // at 1, as in the following:
00386 //
00387 // <srcblock> 
00388 //      REAL    X(n)
00389 //      COMPLEX Y(n/2 + 1)
00390 // </srcblock> 
00391 //
00392 // <example>
00393 // These examples use the table and workspace sizes appropriate to Origin
00394 // series.
00395 //
00396 // Example 1:  Initialize the complex array TABLE in preparation for
00397 // doing an FFT of size 1024.  In this case only the arguments isign, n,
00398 // and table are used. You can use dummy arguments or zeros for the other
00399 // arguments in the subroutine call.
00400 //
00401 // <srcblock> 
00402 //      REAL TABLE(15 + 1024)
00403 //      CALL SCFFT(0, 1024, 0.0, DUMMY, DUMMY, TABLE, DUMMY, 0)
00404 // </srcblock> 
00405 //
00406 // Example 2:  X is a real array of dimension (0:1023), and Y is a
00407 // complex array of dimension (0:512).  Take the FFT of X and store the
00408 // results in Y.  Before taking the FFT, initialize the TABLE array, as
00409 // in example 1.
00410 //
00411 // <srcblock> 
00412 //      REAL X(0:1023)
00413 //      COMPLEX Y(0:512)
00414 //      REAL TABLE(15 + 1024)
00415 //      REAL WORK(1024)
00416 //      CALL SCFFT(0, 1024, 1.0, X, Y, TABLE, WORK, 0)
00417 //      CALL SCFFT(1, 1024, 1.0, X, Y, TABLE, WORK, 0)
00418 // </srcblock> 
00419 //
00420 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
00421 // and store it back in X.  The scale factor 1/1024 is used.  Assume that
00422 // the TABLE array is initialized already.
00423 //
00424 // <srcblock> 
00425 //      CALL CSFFT(-1, 1024, 1.0/1024.0, Y, X, TABLE, WORK, 0)
00426 // </srcblock> 
00427 //
00428 // Example 4:  Perform the same computation as in example 2, but assume
00429 // that the lower bound of each array is 1, rather than 0.  The
00430 // subroutine calls are not changed.
00431 //
00432 // <srcblock> 
00433 //      REAL X(1024)
00434 //      COMPLEX Y(513)
00435 //      CALL SCFFT(0, 1024, 1.0, X, Y, TABLE, WORK, 0)
00436 //      CALL SCFFT(1, 1024, 1.0, X, Y, TABLE, WORK, 0)
00437 // </srcblock> 
00438 //
00439 // Example 5:  Perform the same computation as in example 4, but
00440 // equivalence the input and output arrays to save storage space.  Assume
00441 // that the TABLE array is initialized already.
00442 //
00443 // <srcblock> 
00444 //      REAL X(1024)
00445 //      COMPLEX Y(513)
00446 //      EQUIVALENCE ( X(1), Y(1) )
00447 //      CALL SCFFT(1, 1024, 1.0, X, Y, TABLE, WORK, 0)
00448 // </srcblock> 
00449 // </example>
00450 //
00451 // Input parameters:
00452 // <dl compact>
00453 // <dt><b>isign</b>
00454 // <dd>    Integer.
00455 //         Specifies whether to initialize the table array or to do the
00456 //         forward or inverse Fourier transform, as follows:
00457 //
00458 //         If isign = 0, the routine initializes the table array and
00459 //         returns.  In this case, the only arguments used or checked
00460 //         are isign, n, and table.
00461 //
00462 //         If isign = +1 or -1, the value of isign is the sign of the
00463 //         exponent used in the FFT formula.
00464 // <dt><b>n</b>
00465 // <dd>    Integer.
00466 //         Size of transform.  If n <= 2, <src>scfft/dzfft</src>
00467 //         returns without calculating the transform.
00468 // <dt><b>scale</b>
00469 // <dd>    Scale factor.
00470 //         <src>scfft</src>: real. 
00471 //         <src>dzfft</src>: double precision.
00472 //         <src>csfft</src>: real.
00473 //         <src>zdfft</src>: double precision.
00474 //         Each element of the output array is multiplied by scale
00475 //         after taking the Fourier transform, as defined by the previous
00476 //         formula.
00477 // <dt><b>x</b>
00478 // <dd>    Input array of values to be transformed.
00479 //         <src>scfft</src>: real array of dimension (0:n-1).
00480 //         <src>dzfft</src>: double precision array of dimension (0:n-1).
00481 //         <src>csfft</src>: complex array of dimension (0:n/2).
00482 //         <src>zdfft</src>: double complex array of dimension (0:n/2).
00483 // <dt><b>isys</b>
00484 // <dd>    Integer array of dimension (0:isys(0)).
00485 //         Use isys to specify certain processor-specific parameters or
00486 //         options.  The first element of the array specifies how many
00487 //         more elements are in the array.
00488 //
00489 //         If isys(0) = 0, the default values of such parameters are
00490 //         used.  In this case, you can specify the argument value as
00491 //         the scalar integer constant 0.  If isys(0) > 0, isys(0)
00492 //         gives the upper bound of the isys array; that is, if
00493 //         il = isys(0), user-specified parameters are expected in
00494 //         isys(1) through isys(il).
00495 // </dl>
00496 // Output parameters:
00497 // <dl compact>
00498 // <dt><b>y</b>
00499 // <dd>    Output array of transformed values.
00500 //         <src>scfft</src>: complex array of dimension (0:n/2).
00501 //         <src>dzfft</src>: double complex array of dimension (0:n/2).
00502 //         <src>csfft</src>: real array of dimension (0:n-1).
00503 //         <src>zdfft</src>: double precision array of dimension (0:n-1).
00504 //
00505 //         The output array, y, is the FFT of the the input array, x,
00506 //         computed according to the preceding formula.  The output
00507 //         array may be equivalenced to the input array in the calling
00508 //         program.  Be careful when dimensioning the arrays, in this
00509 //         case, to allow for the fact that the complex array contains
00510 //         two (real) words more than the real array.
00511 // <dt><b>table</b>
00512 // <dd>    Real array; dimension n+15.
00513 //
00514 //         Table of factors and trigonometric functions.
00515 //
00516 //         If isign = 0, the table array is initialized to contain
00517 //         trigonometric tables needed to compute an FFT of size n.
00518 //
00519 //         If isign = +1 or -1, the values in table are assumed to be
00520 //         initialized already by a prior call with isign = 0.
00521 // <dt><b>work</b>
00522 // <dd>    Real array; dimension n.
00523 //
00524 //         Work array used for intermediate calculations.  Its address
00525 //         space must be different from that of the input and output
00526 //         arrays.
00527 // </dl>
00528 // <group>
00529 static void scfft(Int isign, Int n, Float scale, Float* x,
00530                   Complex* y, Float* table, Float* work, Int isys); 
00531 static void scfft(Int isign, Int n, Double scale, Double* x,
00532                   DComplex* y, Double* table, Double* work, Int isys);
00533 static void dzfft(Int isign, Int n, Double scale, Double* x,
00534                   DComplex* y, Double* table, Double* work, Int isys);
00535 static void csfft(Int isign, Int n, Float scale, Complex* x,
00536                   Float* y, Float* table, Float* work, Int isys); 
00537 static void csfft(Int isign, Int n, Double scale, DComplex* x, 
00538                   Double* y, Double* table, Double* work, Int isys);
00539 static void zdfft(Int isign, Int n, Double scale, DComplex* x, 
00540                   Double* y, Double* table, Double* work, Int isys);
00541 // </group>
00542 
00543 // <src>ccfftm/zzfftm</src> computes the FFT of each column of the
00544 // complex matrix x, and stores the results in the columns of complex
00545 // matrix y.
00546 //
00547 // Suppose the arrays are dimensioned as follows:
00548 //
00549 // <srcblock>
00550 //      COMPLEX X(0:ldx-1, 0:lot-1)
00551 //      COMPLEX Y(0:ldy-1, 0:lot-1)
00552 //
00553 // where ldx >= n, ldy >= n.
00554 // </srcblock>
00555 //
00556 // Then column L of the output array is the FFT of column L of the
00557 // input array, using the following formula for the FFT:
00558 //
00559 // <srcblock>
00560 //                        n-1
00561 //      Y(k, L) = scale * Sum [ X(j)*w**(isign*j*k) ]
00562 //                        j=0
00563 //      for k = 0, ..., n-1
00564 //          L = 0, ..., lot-1
00565 //      where:
00566 //          w = exp(2*pi*i/n),
00567 //          i = + sqrt(-1),
00568 //          pi = 3.14159...,
00569 //          isign = +1 or -1
00570 //          lot = the number of columns to transform
00571 // </srcblock>
00572 //
00573 // Different authors use different conventions for which of the
00574 // transforms, isign = +1 or isign = -1, is the forward or inverse
00575 // transform, and what the scale factor should be in either case.  You
00576 // can make this routine compute any of the various possible definitions,
00577 // however, by choosing the appropriate values for isign and scale.
00578 //
00579 // The relevant fact from FFT theory is this:  If you take the FFT with
00580 // any particular values of isign and scale, the mathematical inverse
00581 // function is computed by taking the FFT with -isign and 1/(n * scale).
00582 // In particular, if you use isign = +1 and scale = 1.0 for the forward
00583 // FFT, you can compute the inverse FFT by using the following:  isign =
00584 // -1 and scale = 1.0/n.
00585 //
00586 // This section contains information about the algorithm for these
00587 // routines, the initialization of the table array, the declaration of
00588 // dimensions for x and y arrays, some performance tips, and some
00589 // implementation-dependent details.
00590 //
00591 // <h3>Algorithm</h3>
00592 // These routines use decimation-in-frequency type FFT.  It takes the FFT
00593 // of the columns and vectorizes the operations along the rows of the
00594 // matrix.  Thus, the vector length in the calculations depends on the
00595 // row size, and the strides for vector loads and stores are the leading
00596 // dimensions, ldx and ldy.
00597 //
00598 // <h3>Initialization</h3>
00599 // The table array stores the trigonometric tables used in calculation of
00600 // the FFT.  You must initialize the table array by calling the routine
00601 // with isign = 0 prior to doing the transforms.  If the value of the
00602 // problem size, n, does not change, table does not have to be
00603 // reinitialized.
00604 //
00605 // <h3>Dimensions</h3>
00606 // In the preceding description, it is assumed that array subscripts were
00607 // zero-based, as is customary in FFT applications.  Thus, the input and
00608 // output arrays are declared as follows:
00609 //
00610 // <srcblock> 
00611 //      COMPLEX X(0:ldx-1, 0:lot-1)
00612 //      COMPLEX Y(0:ldy-1, 0:lot-1)
00613 // </srcblock> 
00614 //
00615 // The calling sequence does not have to change, however, if you prefer
00616 // to use the more customary Fortran style with subscripts starting at 1.
00617 // The same values of ldx and ldy would be passed to the subroutine even
00618 // if the input and output arrays were dimensioned as follows:
00619 //
00620 // <srcblock> 
00621 //      COMPLEX X(ldx, lot)
00622 //      COMPLEX Y(ldy, lot)
00623 // </srcblock> 
00624 //
00625 // <example>
00626 // Example 1:  Initialize the TABLE array in preparation for doing an FFT
00627 // of size 128.  Only the isign, n, and table arguments are used in this
00628 // case.  You can use dummy arguments or zeros for the other arguments in
00629 // the subroutine call.
00630 //
00631 // <srcblock> 
00632 //       REAL TABLE(30 + 256)
00633 //       CALL CCFFTM(0, 128, 0, 0., DUMMY, 1, DUMMY, 1, TABLE, DUMMY, 0)
00634 // </srcblock> 
00635 //
00636 // Example 2:  X and Y are complex arrays of dimension (0:128) by (0:55).
00637 // The first 128 elements of each column contain data.  For performance
00638 // reasons, the extra element forces the leading dimension to be an odd
00639 // number.  Take the FFT of the first 50 columns of X and store the
00640 // results in the first 50 columns of Y.  Before taking the FFT,
00641 // initialize the TABLE array, as in example 1.
00642 //
00643 // <srcblock> 
00644 //       COMPLEX X(0:128, 0:55)
00645 //       COMPLEX Y(0:128, 0:55)
00646 //       REAL TABLE(30 + 256)
00647 //       REAL WORK(256)
00648 //       ...
00649 //       CALL CCFFTM(0, 128, 50, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
00650 //       CALL CCFFTM(1, 128, 50, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
00651 // </srcblock> 
00652 //
00653 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
00654 // and store it back in X.  The scale factor 1/128 is used.  Assume that
00655 // the TABLE array is already initialized.
00656 //
00657 // <srcblock> 
00658 //       CALL CCFFTM(-1, 128, 50, 1./128., Y, 129, X, 129, TABLE,WORK,0)
00659 // </srcblock> 
00660 //
00661 // Example 4:  Perform the same computation as in example 2, but assume
00662 // that the lower bound of each array is 1, rather than 0.  The
00663 // subroutine calls are not changed.
00664 //
00665 // <srcblock> 
00666 //       COMPLEX X(129, 55)
00667 //       COMPLEX Y(129, 55)
00668 //       ...
00669 //       CALL CCFFTM(0, 128, 50, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
00670 //       CALL CCFFTM(1, 128, 50, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
00671 // </srcblock> 
00672 //
00673 // Example 5:  Perform the same computation as in example 4, but put the
00674 // output back in array X to save storage space.  Assume that the TABLE
00675 // array is already initialized.
00676 //
00677 // <srcblock> 
00678 //       COMPLEX X(129, 55)
00679 //       ...
00680 //       CALL CCFFTM(1, 128, 50, 1.0, X, 129, X, 129, TABLE, WORK, 0)
00681 // </srcblock> 
00682 //
00683 // </example>
00684 //
00685 // Input parameters:
00686 // <dl compact>
00687 // <dt><b>isign</b>
00688 // <dd>    Integer.
00689 //         Specifies whether to initialize the table array or to do the
00690 //         forward or inverse Fourier transform, as follows:
00691 //
00692 //         If isign = 0, the routine initializes the table array and
00693 //         returns.  In this case, the only arguments used or checked
00694 //         are isign, n, and table.
00695 //
00696 //         If isign = +1 or -1, the value of isign is the sign of the
00697 //         exponent used in the FFT formula.
00698 // <dt><b>n</b>
00699 // <dd>    Integer.  
00700 //         Size of each transform (the number of elements in each
00701 //         column of the input and output matrix to be transformed).
00702 //         Performance depends on the value of n, as explained above.  
00703 //         n >= 0; if n = 0, the routine returns.
00704 // <dt><b>lot</b>
00705 // <dd>    Integer.
00706 //         The number of transforms to be computed (lot size).  This is
00707 //         the number of elements in each row of the input and output
00708 //         matrix.  lot >= 0.  If lot = 0, the routine returns.
00709 // <dt><b>scale</b>
00710 // <dd>    Scale factor.
00711 //         <src>ccfftm</src>: real. 
00712 //         <src>zzfftm</src>: double precision.
00713 //         Each element of the output array is multiplied by scale
00714 //         factor after taking the Fourier transform, as defined
00715 //         previously.
00716 // <dt><b>x</b>
00717 // <dd>    Array of dimension (0:ldx-1, 0:n2-1).
00718 //         <src>ccfftm</src>: real array.
00719 //         <src>zzfftm</src>: double precision array.
00720 //         Input array of values to be transformed.
00721 // <dt><b>ldx</b>
00722 // <dd>    The number of rows in the x array, as it was declared in the
00723 //         calling program (the leading dimension of X).  ldx >= MAX(n, 1).
00724 // <dt><b>ldy</b>
00725 // <dd>    Integer.
00726 //         The number of rows in the y array, as it was declared in the
00727 //         calling program (the leading dimension of y).  ldy >= MAX(n,
00728 //         1).
00729 // <dt><b>isys</b>
00730 // <dd>    Integer array of dimension (0:isys(0)).
00731 //         The first element of the array specifies how many more
00732 //         elements are in the array.  Use isys to specify certain
00733 //         processor-specific parameters or options.
00734 //
00735 //         If isys(0) = 0, the default values of such parameters are
00736 //         used.  In this case, you can specify the argument value as
00737 //         the scalar integer constant 0.
00738 //
00739 //         If isys(0) > 0, isys(0) gives the upper bound of the isys
00740 //         array.  Therefore, if il = isys(0), user-specified
00741 //         parameters are expected in isys(1) through isys(il).
00742 // </dl>
00743 // Output parameters:
00744 // <dl compact>
00745 // <dt><b>y</b>
00746 // <dd>    Array of dimension (0:ldy-1, 0:lot-1).
00747 //         <src>ccfftm</src>: complex array.
00748 //         <src>zzfftm</src>: double complex array.
00749 //         Output array of transformed values.  Each column of the
00750 //         output array, y, is the FFT of the corresponding column of
00751 //         the input array, x, computed according to the preceding
00752 //         formula.
00753 //
00754 //         The output array may be the same as the input array. In that
00755 //         case, the transform is done in place.  The input array is
00756 //         overwritten with the transformed values.  In this case, it
00757 //         is necessary that ldx = ldy.
00758 // <dt><b>table</b>
00759 // <dd>    Real array; dimension (30 + 2n).
00760 //         Table of factors and trigonometric functions.
00761 //
00762 //         If isign = 0, the routine initializes table (table is output
00763 //         only).
00764 //
00765 //         If isign = +1 or -1, the values in table are assumed to be
00766 //         initialized already by a prior call with isign = 0 (table is
00767 //         input only).
00768 // <dt><b>work</b>
00769 // <dd>    Real array; dimension 2n.
00770 //         Work array.  This is a scratch array used for intermediate
00771 //         calculations.  Its address space must be different from that
00772 //         of the input and output arrays.
00773 // </dl>
00774 // <group>
00775 static void ccfftm(Int isign, Int n, Int lot, Float scale, Complex*
00776                    x, Int ldx, Complex* y, Int ldy, Float* table,
00777                    Float* work, Int isys); 
00778 static void zzfftm(Int isign, Int n, Int lot, Double scale, DComplex*
00779                    x, Int ldx, DComplex* y, Int ldy, Double* table,
00780                    Double* work, Int isys);
00781 // </group>
00782 
00783 // <src>scfftm/dzfftm</src> computes the FFT of each column of the real matrix
00784 // X, and it stores the results in the corresponding column of the complex
00785 // matrix Y.  <src>csfftm/zdfftm</src> computes the corresponding inverse
00786 // transforms.
00787 //
00788 // In FFT applications, it is customary to use zero-based subscripts; the
00789 // formulas are simpler that way.  First, the function of <src>scfftm</src> is
00790 // described.  Suppose that the arrays are dimensioned as follows:
00791 //
00792 // <srcblock>
00793 //      REAL    X(0:ldx-1, 0:lot-1)
00794 //      COMPLEX Y(0:ldy-1, 0:lot-1)
00795 //
00796 // where ldx >= n, ldy >= n/2 + 1.
00797 // </srcblock>
00798 //
00799 // Then column L of the output array is the FFT of column L of the input
00800 // array, using the following formula for the FFT:
00801 //
00802 // <srcblock>
00803 //                    n-1
00804 // Y(k, L) = scale *  Sum  [ X(j, L)*w**(isign*j*k) ]
00805 //                    j=0
00806 //
00807 // for k = 0, ..., n/2
00808 //     L = 0, ..., lot-1 where:
00809 //     w = exp(2*pi*i/n),
00810 //     i = + sqrt(-1)
00811 //     pi = 3.14159...,
00812 //     isign = +1 or -1,
00813 //     lot = the number of columns to transform
00814 // </srcblock>
00815 //
00816 // Different authors use different conventions for which transform
00817 // (isign = +1 or isign = -1) is used in the real-to-complex case, and
00818 // what the scale factor should be.  Some adopt the convention that isign
00819 // = 1 for the real-to-complex transform, and isign = -1 for the
00820 // complex-to-real inverse.  Others use the opposite convention.  You can
00821 // make these routines compute any of the various possible definitions,
00822 // however, by choosing the appropriate values for isign and scale.
00823 //
00824 // The relevant fact from FFT theory is this:  If you use <src>scfftm</src> to
00825 // take the real-to-complex FFT, using any particular values of isign and
00826 // scale, the mathematical inverse function is computed by using 
00827 // <src>csfftm</src> with -isign and 1/ (n*scale).  In particular, if you call
00828 // <src>scfftm</src> with isign = +1 and scale = 1.0, you can use
00829 // <src>csfftm</src> to compute the inverse complex-to-real FFT by using isign
00830 // = -1 and scale = 1.0/n.
00831 //
00832 // <h3>Real-to-complex FFTs</h3>
00833 // Notice in the preceding formula that there are n real input values and
00834 // (n/2) + 1 complex output values for each column.  This property is
00835 // characteristic of real-to-complex FFTs.
00836 //
00837 // The mathematical definition of the Fourier transform takes a sequence
00838 // of n complex values and transforms it to another sequence of n complex
00839 // values.  A complex-to-complex FFT routine, such as <src>ccfftm</src>, will
00840 // take n complex input values and produce n complex output values.  In fact,
00841 // one easy way to compute a real-to-complex FFT is to store the input
00842 // data x in a complex array, then call routine <src>ccfftm</src> to compute
00843 // the FFT.  You get the same answer when using the <src>scfftm</src> routine.
00844 //
00845 // A separate real-to-complex FFT routine is more efficient than the
00846 // equivalent complex-to-complex routine.  Because the input data is
00847 // real, you can make use of this fact to save almost half of the
00848 // computational work.  According to the theory of Fourier transforms,
00849 // for real input data, you have to compute only the first n/2 + 1
00850 // complex output values in each column, because the second half of the
00851 // FFT values in each column can be computed from the first half of the
00852 // values by the simple formula:
00853 //
00854 // <srcblock>
00855 //      Y    = conjgY          for n/2 <= k <= n-1
00856 //       k,L          n-k, L
00857 //
00858 // where the notation conjg(z) represents the complex conjugate of z.
00859 // </srcblock>
00860 //
00861 // In fact, in many applications, the second half of the complex output
00862 // data is never explicitly computed or stored.  Likewise, you must
00863 // supply only the first half of the complex data in each column has to
00864 // be supplied for the complex-to-real FFT.
00865 //
00866 // Another implication of FFT theory is that for real input data, the
00867 // first output value in each column, Y(0, L), will always be a real
00868 // number; therefore, the imaginary part will always be 0.  If n is an
00869 // even number, Y(n/2, L) will also be real and have 0 imaginary parts.
00870 //
00871 // <h3>Complex-to-real FFTs</h3>
00872 // Consider the complex-to-real case.  The effect of the computation is
00873 // given by the preceding formula, but with X complex and Y real.
00874 //
00875 // In general, the FFT transforms a complex sequence into a complex
00876 // sequence; however, in a certain application you may know the output
00877 // sequence is real, perhaps because the complex input sequence was the
00878 // transform of a real sequence.  In this case, you can save about half
00879 // of the computational work.
00880 //
00881 // According to the theory of Fourier transforms, for the output
00882 // sequence, Y, to be a real sequence, the following identity on the
00883 // input sequence, X, must be true:
00884 //
00885 // <srcblock>
00886 //      X    = conjgX         for n/2 <= k <= n-1
00887 //       k,L          n-k,L
00888 // And, in fact, the following input values
00889 //
00890 //      X    for k > n/2
00891 //       k,L
00892 // do not have to be supplied, because they can be inferred from the
00893 // first half of the input.
00894 // </srcblock>
00895 //
00896 // Thus, in the complex-to-real routine, CSFFTM, the arrays can be
00897 // dimensioned as follows:
00898 //
00899 // <srcblock>
00900 //      COMPLEX X(0:ldx-1, 0:lot-1)
00901 //      REAL    Y(0:ldy-1, 0:lot-1)
00902 //
00903 // where ldx >= n/2 + 1, ldy >= n.
00904 // </srcblock>
00905 //
00906 // In each column, there are (n/2) + 1 complex input values and n real
00907 // output values.  Even though only (n/2) + 1 input values are supplied,
00908 // the size of the transform is still n in this case, because implicitly
00909 // the FFT formula for a sequence of length n is used.
00910 //
00911 // Another implication of the theory is that X(0, L) must be a real
00912 // number (that is, must have zero imaginary part).  If n is an even
00913 // number, X(n/2, L) must also be real.  Routine CSFFTM assumes that each
00914 // of these values is real; if a nonzero imaginary part is given, it is
00915 // ignored.
00916 //
00917 // <h3>Table Initialization</h3>
00918 // The table array contains the trigonometric tables used in calculation
00919 // of the FFT.  You must initialize this table by calling the routine
00920 // with isign = 0 prior to doing the transforms.  table does not have to
00921 // be reinitialized if the value of the problem size, n, does not change.
00922 //
00923 // <h3>Dimensions</h3>
00924 // In the preceding description, it is assumed that array subscripts were
00925 // zero-based, as is customary in FFT applications.  Thus, the input and
00926 // output arrays are declared (for SCFFTM):
00927 //
00928 // <srcblock>
00929 //      REAL    X(0:ldx-1, 0:lot-1)
00930 //      COMPLEX Y(0:ldy-1, 0:lot-1)
00931 // </srcblock>
00932 //
00933 // No change is made in the calling sequence, however, if you prefer to
00934 // use the more customary Fortran style with subscripts starting at 1.
00935 // The same values of ldx and ldy would be passed to the subroutine even
00936 // if the input and output arrays were dimensioned as follows:
00937 //
00938 // <srcblock>
00939 //      REAL    X(ldx, lot)
00940 //      COMPLEX Y(ldy, lot)
00941 // </srcblock>
00942 
00943 // </example>
00944 // Example 1:  Initialize the complex array TABLE in preparation for
00945 // doing an FFT of size 128.  In this case only the isign, n, and table
00946 // arguments are used; you may use dummy arguments or zeros for the other
00947 // arguments in the subroutine call.
00948 //
00949 // <srcblock>
00950 //       REAL TABLE(15 + 128)
00951 //       CALL SCFFTM(0, 128, 1, 0.0, DUMMY, 1, DUMMY, 1,
00952 //      &  TABLE, DUMMY, 0)
00953 // </srcblock>
00954 //
00955 // Example 2:  X is a real array of dimension (0:128, 0:55), and Y is a
00956 // complex array of dimension (0:64, 0:55).  The first 128 elements in
00957 // each column of X contain data; the extra element forces an odd leading
00958 // dimension.  Take the FFT of the first 50 columns of X and store the
00959 // results in the first 50 columns of Y.  Before taking the FFT,
00960 // initialize the TABLE array, as in example 1.
00961 //
00962 // <srcblock>
00963 //       REAL    X(0:128, 0:55)
00964 //       COMPLEX Y(0:64,  0:55)
00965 //       REAL    TABLE(15 + 128)
00966 //       REAL    WORK((128)
00967 //       ...
00968 //       CALL SCFFTM(0, 128, 50, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
00969 //       CALL SCFFTM(1, 128, 50, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
00970 // </srcblock>
00971 //
00972 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
00973 // and store it back in X.  The scale factor 1/128 is used.  Assume that
00974 // the TABLE array is initialized already.
00975 //
00976 // <srcblock>
00977 //       CALL CSFFTM(-1, 128, 50, 1.0/128.0, Y, 65, X, 129,
00978 //      &  TABLE, WORK, 0)
00979 // </srcblock>
00980 //
00981 // Example 4:  Perform the same computation as in example 2, but assume
00982 // that the lower bound of each array is 1, rather than 0.  No change is
00983 // made in the subroutine calls.
00984 //
00985 // <srcblock>
00986 //       REAL    X(129, 56)
00987 //       COMPLEX Y(65, 56)
00988 //       ...
00989 //       CALL SCFFTM(0, 128, 50, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
00990 //       CALL SCFFTM(1, 128, 50, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
00991 // </srcblock>
00992 //
00993 // Example 5:  Perform the same computation as in example 4, but
00994 // equivalence the input and output arrays to save storage space.  In
00995 // this case, a row must be added to X, because it is equivalenced to a
00996 // complex array.  The leading dimension of X is two times an odd number;
00997 // therefore, memory bank conflicts are minimal.  Assume that TABLE is
00998 // initialized already.
00999 //
01000 // <srcblock>
01001 //       REAL    X(130, 56)
01002 //       COMPLEX Y(65, 56)
01003 //       EQUIVALENCE ( X(1, 1), Y(1, 1) )
01004 //       ...
01005 //       CALL SCFFTM(1, 128, 50, 1.0, X, 130, Y, 65, TABLE, WORK, 0)
01006 // </srcblock>
01007 // </example>
01008 //
01009 // Input parameters:
01010 // <dl compact>
01011 // <dt><b>isign</b>
01012 // <dd>    Integer.
01013 //         Specifies whether to initialize the table array or to do the
01014 //         forward or inverse Fourier transform, as follows:
01015 //
01016 //         If isign = 0, the routine initializes the table array and
01017 //         returns.  In this case, the only arguments used or checked
01018 //         are isign, n, and table.
01019 //
01020 //         If isign = +1 or -1, the value of isign is the sign of the
01021 //         exponent used in the FFT formula.
01022 // <dt><b>n</b>
01023 // <dd>    Integer.
01024 //         Size of the transforms (the number of elements in each
01025 //         column of the input and output matrix to be transformed).
01026 //         If n is not positive, <src>scfftm</src> or <src>csfftm</src> returns
01027 //         without computing a transforms.
01028 // <dt><b>lot</b>
01029 // <dd>    Integer.
01030 //         The number of transforms to be computed (or "lot size").
01031 //         This is the number of elements in each row of the input and
01032 //         output matrix.  If lot is not positive, <src>csfftm</src> or 
01033 //         <src>scfftm</src> returns without computing a transforms.
01034 // <dt><b>scale</b>
01035 // <dd>    Scale factor.
01036 //         <src>scfftm</src>: real. 
01037 //         <src>dzfftm</src>: double precision.
01038 //         <src>csfftm</src>: real.
01039 //         <src>zdfftm</src>: double precision.
01040 //         Each element of the output array is multiplied by scale
01041 //         after taking the transform, as defined in the preceding
01042 //         formula.
01043 // <dt><b>x</b>
01044 // <dd>    Input array of values to be transformed. Dimension (0:ldx-1,
01045 //         0:lot-1).
01046 //         <src>scfftm</src>: real array.
01047 //         <src>dzfftm</src>: double precision array.
01048 //         <src>csfftm</src>: complex array.
01049 //         <src>zdfftm</src>: double complex array.
01050 // <dt><b>ldx</b>
01051 // <dd>    Integer.
01052 //         The number of rows in the x array, as it was declared in the
01053 //         calling program.  That is, the leading dimension of x.
01054 //         <src>scfftm, dzfftm</src>:  ldx >= MAX(n, 1).
01055 //         <src>csfftm, zdfftm</src>:  ldx >= MAX(n/2 + 1, 1).
01056 // <dt><b>ldy</b>
01057 // <dd>    Integer.
01058 //         The number of rows in the y array, as it was declared in the
01059 //         calling program (the leading dimension of y).
01060 //         <src>scfftm, dzfftm</src>:  ldy >= MAX(n/2 + 1, 1).
01061 //         <src>csfftm, zdfftm</src>:  ldy >= MAX(n, 1).
01062 // <dt><b>isys</b>
01063 // <dd>    Integer array of dimension (0:isys(0)).
01064 //         The first element of the array specifies how many more
01065 //         elements are in the array.  Use isys to specify certain
01066 //         processor-specific parameters or options.
01067 //
01068 //         If isys(0) = 0, the default values of such parameters are
01069 //         used.  In this case, you can specify the argument value as
01070 //         the scalar integer constant 0.
01071 //
01072 //         If isys(0) > 0, isys(0) gives the upper bound of the isys
01073 //         array.  Therefore, if il = isys(0), user-specified
01074 //         parameters are expected in isys(1) through isys(il).
01075 // </dl>
01076 // Output parameters:
01077 // <dl compact>
01078 // <dt><b>y</b>
01079 // <dd>    Output array of transformed values.  Dimension (0:ldy-1,
01080 //         0:lot-1).
01081 //         <src>scfftm</src>: complex array.
01082 //         <src>dzfftm</src>: double complex array.
01083 //         <src>csfftm</src>: real array.
01084 //         <src>zdfftm</src>: double precision array.
01085 //
01086 //         Each column of the output array, y, is the FFT of the
01087 //         corresponding column of the input array, x, computed
01088 //         according to the preceding formula.  The output array may be
01089 //         equivalenced to the input array. In that case, the transform
01090 //         is done in place and the input array is overwritten with the
01091 //         transformed values.  In this case, the following conditions
01092 //         on the leading dimensions must hold:
01093 //
01094 //         <src>scfftm, dzfftm</src>:  ldx = 2ldy.
01095 //         <src>csfftm, zdfftm</src>:  ldy = 2ldx.
01096 // <dt><b>table</b>
01097 // <dd>    Real array; dimension (15 + n).
01098 //         Table of factors and trigonometric functions.
01099 //         This array must be initialized by a call to <src>scfftm</src> (or
01100 //         <src>csfftm</src>) with isign = 0.
01101 //
01102 //         If isign = 0, table is initialized to contain trigonometric
01103 //         tables needed to compute an FFT of length n.
01104 // <dt><b>work</b>
01105 // <dd>    Real array; dimension n.
01106 //         Work array used for intermediate calculations.  Its address
01107 //         space must be different from that of the input and output
01108 //         arrays.
01109 // </dl>
01110 // <group>
01111 static void scfftm(Int isign, Int n, Int lot, Float scale, Float*
01112                    x, Int ldx, Complex* y, Int ldy, Float* table,
01113                    Float* work, Int isys); 
01114 static void dzfftm(Int isign, Int n, Int lot, Double scale, Double*
01115                    x, Int ldx, DComplex* y, Int ldy, Double* table,
01116                    Double* work, Int isys); 
01117 static void csfftm(Int isign, Int n, Int lot, Float scale, Complex*
01118                    x, Int ldx, Float* y, Int ldy, Float* table,
01119                    Float* work, Int isys); 
01120 static void zdfftm(Int isign, Int n, Int lot, Double scale, DComplex*
01121                    x, Int ldx, Double* y, Int ldy, Double* table,
01122                    Double* work, Int isys);
01123 // </group>
01124 
01125 // These routines compute the two-dimensional complex Fast Fourier
01126 // Transform (FFT) of the complex matrix x, and store the results in the
01127 // complex matrix y.  <src>ccfft2d</src> does the complex-to-complex
01128 // transform and <src>zzfft</src> does the same for double
01129 // precision arrays.
01130 //
01131 // In FFT applications, it is customary to use zero-based subscripts; the
01132 // formulas are simpler that way.  Suppose that the arrays are
01133 // dimensioned as follows:
01134 //
01135 // <srcblock> 
01136 //      COMPLEX X(0:n1-1, 0:n2-1)
01137 //      COMPLEX Y(0:n1-1, 0:n2-1)
01138 // </srcblock> 
01139 //
01140 // These routines compute the formula:
01141 //
01142 // <srcblock> 
01143 //                        n2-1  n1-1
01144 //    Y(k1, k2) = scale * Sum   Sum [ X(j1, j2)*w1**(j1*k1)*w2**(j2*k2) ]
01145 //                        j2=0  j1=0
01146 //
01147 //    for k1 = 0, ..., n1-1
01148 //        k2 = 0, ..., n2-1
01149 //
01150 //    where:
01151 //        w1 = exp(isign*2*pi*i/n1)
01152 //        w2 = exp(isign*2*pi*i/n2)
01153 //        i = + sqrt(-1)
01154 //        pi = 3.14159...,
01155 //        isign = +1 or -1
01156 // </srcblock> 
01157 //
01158 // Different authors use different conventions for which of the
01159 // transforms, isign = +1 or isign = -1, is the forward or inverse
01160 // transform, and what the scale factor should be in either case.  You
01161 // can make this routine compute any of the various possible definitions,
01162 // however, by choosing the appropriate values for isign and scale.
01163 //
01164 // The relevant fact from FFT theory is this:  If you take the FFT with
01165 // any particular values of isign and scale, the mathematical inverse
01166 // function is computed by taking the FFT with -isign and
01167 // 1/(n1*n2*scale).  In particular, if you use isign = +1 and scale = 1.0
01168 // for the forward FFT, you can compute the inverse FFT by using isign =
01169 // -1 and scale = 1.0/(n1*n2).
01170 //
01171 // <h3>Algorithm</h3>
01172 // These routines use a routine very much like <src>ccfftm/zzfftm</src> to do
01173 // multiple FFTs first on all columns in an input matrix and then on all
01174 // of the rows.
01175 //
01176 // <h3>Initialization</h3>
01177 // The table array stores factors of n1 and n2 and also trigonometric
01178 // tables that are used in calculation of the FFT.  This table must be
01179 // initialized by calling the routine with isign = 0.  If the values of
01180 // the problem sizes, n1 and n2, do not change, the table does not have
01181 // to be reinitialized.
01182 //
01183 // <h3>Dimensions</h3>
01184 // In the preceding description, it is assumed that array subscripts were
01185 // zero-based, as is customary in FFT applications.  Thus, the input and
01186 // output arrays are declared as follows:
01187 //
01188 // <srcblock> 
01189 //      COMPLEX X(0:ldx-1, 0:n2-1)
01190 //      COMPLEX Y(0:ldy-1, 0:n2-1)
01191 // </srcblock> 
01192 //
01193 // However, the calling sequence does not change if you prefer to use the
01194 // more customary Fortran style with subscripts starting at 1.  The same
01195 // values of ldx and ldy would be passed to the subroutine even if the
01196 // input and output arrays were dimensioned as follows:
01197 //
01198 // <srcblock> 
01199 //      COMPLEX X(ldx, n2)
01200 //      COMPLEX Y(ldy, n2)
01201 // </srcblock> 
01202 //
01203 // <example>
01204 // All examples here are for Origin series only.
01205 //
01206 // Example 1:  Initialize the TABLE array in preparation for doing a
01207 // two-dimensional FFT of size 128 by 256.  In this case only the isign,
01208 // n1, n2, and table arguments are used; you can use dummy arguments or
01209 // zeros for other arguments.
01210 //
01211 // <srcblock> 
01212 //        REAL TABLE ((30 + 256) + (30 + 512))
01213 //        CALL CCFFT2D (0, 128, 256, 0.0, DUMMY, 1, DUMMY, 1,
01214 //       &  TABLE, DUMMY, 0)
01215 // </srcblock> 
01216 //
01217 // Example 2:  X and Y are complex arrays of dimension (0:128, 0:255).
01218 // The first 128 elements of each column contain data.  For performance
01219 // reasons, the extra element forces the leading dimension to be an odd
01220 // number.  Take the two-dimensional FFT of X and store it in Y.
01221 // Initialize the TABLE array, as in example 1.
01222 //
01223 // <srcblock> 
01224 //       COMPLEX X(0:128, 0:255)
01225 //       COMPLEX Y(0:128, 0:255)
01226 //       REAL     TABLE((30 + 256) + (30 + 512))
01227 //       REAL     WORK(2*128*256)
01228 //       ...
01229 //       CALL CCFFT2D(0, 128, 256, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
01230 //       CALL CCFFT2D(1, 128, 256, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
01231 // </srcblock> 
01232 //
01233 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
01234 // and store it back in X.  The scale factor 1/(128*256) is used.  Assume
01235 // that the TABLE array is already initialized.
01236 //
01237 // <srcblock> 
01238 //       CALL CCFFT2D(-1, 128, 256, 1.0/(128.0*256.0), Y, 129,
01239 //      &  X, 129, TABLE, WORK, 0)
01240 // </srcblock> 
01241 //
01242 // Example 4:  Perform the same computation as in example 2, but assume
01243 // that the lower bound of each array is 1, rather than 0.  The
01244 // subroutine calls are not changed.
01245 //
01246 // <srcblock> 
01247 //       COMPLEX X(129, 256)
01248 //       COMPLEX Y(129, 256)
01249 //       ...
01250 //       CALL CCFFT2D(0, 128, 256, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
01251 //       CALL CCFFT2D(1, 128, 256, 1.0, X, 129, Y, 129, TABLE, WORK, 0)
01252 // </srcblock> 
01253 //
01254 // Example 5:  Perform the same computation as in example 4, but put the
01255 // output back in array X to save storage space.  Assume that the TABLE
01256 // array is already initialized.
01257 //
01258 // <srcblock> 
01259 //       COMPLEX X(129, 256)
01260 //       ...
01261 //       CALL CCFFT2D(1, 128, 256, 1.0, X, 129, X, 129, TABLE, WORK, 0)
01262 // </srcblock> 
01263 // </example>
01264 //
01265 // Input parameters:
01266 // <dl compact>
01267 // <dt><b>isign</b>
01268 // <dd>    Integer.
01269 //         Specifies whether to initialize the table array or to do the
01270 //         forward or inverse transform as follows:
01271 //
01272 //         If isign = 0, the routine initializes the table array and
01273 //         returns.  In this case, the only arguments used or checked
01274 //         are isign, n1, n2, table.
01275 //
01276 //         If isign = +1 or -1, the value of isign is the sign of the
01277 //         exponent used in the FFT formula.
01278 // <dt><b>n1</b>
01279 // <dd>    Integer.
01280 //         Transform size in the first dimension.  If n1 is not
01281 //         positive, the routine returns without performing a
01282 //         transform.
01283 // <dt><b>n2</b>
01284 // <dd>    Integer.
01285 //         Transform size in the second dimension.  If n2 is not
01286 //         positive, the routine returns without performing a
01287 //         transform.
01288 // <dt><b>scale</b>
01289 // <dd>    Scale factor.
01290 //         ccfft2d: real.
01291 //         zzfft2d: double precision.
01292 //         Each element of the output array is multiplied by scale
01293 //         factor after taking the Fourier transform, as defined
01294 //         previously.
01295 // <dt><b>x</b>
01296 // <dd>    Array of dimension (0:ldx-1, 0:n2-1).
01297 //         ccfft2d: complex array.
01298 //         zzfft2d: double complex array.
01299 //         Input array of values to be transformed.
01300 // <dt><b>ldx</b>
01301 // <dd>    Integer.
01302 //         The number of rows in the x array, as it was declared in the
01303 //         calling program (the leading dimension of x).  ldx >=
01304 //         MAX(n1, 1).
01305 // <dt><b>ldy</b>
01306 // <dd>    Integer.
01307 //
01308 //         The number of rows in the y array, as it was declared in the
01309 //         calling program (the leading dimension of y).  ldy >=
01310 //         MAX(n1, 1).
01311 // <dt><b>isys</b>
01312 // <dd>    Algorithm used; value dependent on hardware system.  Currently, no
01313 //         special options are supported; therefore, you must always specify
01314 //         an isys argument as constant 0.  
01315 // </dl>
01316 // Output parameters:
01317 // <dl compact>
01318 // <dt><b>y</b>
01319 // <dd>    Array of dimension (0:ldy-1, 0:n2-1).
01320 //         ccfft2d: complex array.
01321 //         zzfft2d: double complex array.
01322 //         Output array of transformed values.  The output array may be
01323 //         the same as the input array, in which case, the transform is
01324 //         done in place (the input array is overwritten with the
01325 //         transformed values).  In this case, it is necessary that
01326 //         ldx = ldy.
01327 // <dt><b>table</b>
01328 // <dd>    Real array; dimension (30+ 2 * n1) + (30 + 2 * n2).
01329 //
01330 //         Table of factors and trigonometric functions.
01331 //
01332 //         If isign = 0, the routine initializes table (table is output
01333 //         only).
01334 //
01335 //         If isign = +1 or -1, the values in table are assumed to be
01336 //         initialized already by a prior call with isign = 0 (table is
01337 //         input only).
01338 // <dt><b>work</b>
01339 // <dd>    Real array; dimension 2 * (n1*n2).
01340 //
01341 //         Work array.  This is a scratch array used for intermediate
01342 //         calculations.  Its address space must be different from that
01343 //         of the input and output arrays.
01344 // </dl>
01345 // <group>
01346 static void ccfft2d(Int isign, Int n1, Int n2, Float scale, Complex*
01347                     x, Int ldx, Complex* y, Int ldy, Float* table,
01348                     Float* work, Int isys); 
01349 static void zzfft2d(Int isign, Int n1, Int n2, Double scale, DComplex*
01350                     x, Int ldx, DComplex* y, Int ldy, Double* table,
01351                     Double* work, Int isys);
01352 // </group>
01353  
01354 // <src>scfft2d/dzfft2d</src> computes the two-dimensional Fast Fourier
01355 // Transform (FFT) of the real matrix x, and it stores the results in the
01356 // complex matrix y.  <src>csfft2d/zdfft2d</src> computes the corresponding
01357 // inverse transform.
01358 //
01359 // In FFT applications, it is customary to use zero-based subscripts; the
01360 // formulas are simpler that way.  First the function of <src>scfft2d</src> is
01361 // described.  Suppose the arrays are dimensioned as follows:
01362 //
01363 // <srcblock>
01364 //      REAL    X(0:ldx-1, 0:n2-1)
01365 //      COMPLEX Y(0:ldy-1, 0:n2-1)
01366 //
01367 // where ldx >= n1 ldy >= (n1/2) + 1.
01368 // </srcblock>
01369 //
01370 // <src>scfft2d</src> computes the formula:
01371 //
01372 // <srcblock>
01373 //                         n2-1 n1-1
01374 //     Y(k1, k2) = scale * Sum  Sum [ X(j1, j2)*w1**(j1*k1)*w2**(j2*k2) ]
01375 //                         j2=0 j1=0
01376 //
01377 //     for k1 = 0, ..., n1/2 + 1
01378 //         k2 = 0, ..., n2-1
01379 //
01380 //     where:
01381 //         w1 = exp(isign*2*pi*i/n1)
01382 //         w2 = exp(isign*2*pi*i/n2)
01383 //         i  = + sqrt(-1)
01384 //         pi = 3.14159...,
01385 //         isign = +1 or -1
01386 // </srcblock>
01387 //
01388 // Different authors use different conventions for which of the
01389 // transforms, isign = +1 or isign = -1, is the forward or inverse
01390 // transform, and what the scale factor should be in either case.  You
01391 // can make these routines compute any of the various possible
01392 // definitions, however, by choosing the appropriate values for isign and
01393 // scale.
01394 //
01395 // The relevant fact from FFT theory is this:  If you take the FFT with
01396 // any particular values of isign and scale, the mathematical inverse
01397 // function is computed by taking the FFT with -isign and 1/(n1 * n2 *
01398 // scale).  In particular, if you use isign = +1 and scale = 1.0 for the
01399 // forward FFT, you can compute the inverse FFT by using isign = -1 and
01400 // scale = 1.0/(n1 . n2).
01401 //
01402 // <src>scfft2d</src> is very similar in function to <src>ccfft2d</src>, but
01403 // it takes the real-to-complex transform in the first dimension, followed by
01404 // the complex-to-complex transform in the second dimension.
01405 //
01406 // <src>csfft2d</src> does the reverse.  It takes the complex-to-complex FFT
01407 // in the second dimension, followed by the complex-to-real FFT in the first
01408 // dimension.
01409 //
01410 // See the <src>scfft</src> man page for more information about real-to-complex
01411 // and complex-to-real FFTs.  The two-dimensional analog of the conjugate
01412 // formula is as follows:
01413 //
01414 // <srcblock>
01415 //      Y       = conjg Y
01416 //       k , k            n1 - k , n2 - k
01417 //        1   2                 1        2
01418 //
01419 //      for n1/2 <  k  <= n1 - 1
01420 //                   1
01421 //
01422 //         0 <= k  <= n2 - 1
01423 //               2
01424 // where the notation conjg(z) represents the complex conjugate of z.
01425 // </srcblock>
01426 //
01427 // Thus, you have to compute only (slightly more than) half of the output
01428 // values, namely:
01429 //
01430 // <srcblock>
01431 //      Y       for 0 <= k  <= n1/2    0 <= k  <= n2 - 1
01432 //       k , k            1                  2
01433 //        1   2
01434 // </srcblock>
01435 //
01436 // <h3>Algorithm</h3>
01437 // <src>scfft2d</src> uses a routine similar to <src>scfftm</src> to do a
01438 // real-to-complex FFT on the columns, then uses a routine similar to
01439 // <src>ccfftm</src> to do a complex-to-complex FFT on the rows.
01440 //
01441 // <src>csfft2d</src> uses a routine similar to <src>ccfftm</src> to do a
01442 // complex-to-complex FFT on the rows, then uses a routine similar to
01443 // <src>csfftm</src> to do a complex-to-real FFT on the columns.
01444 //
01445 // <h3>Table Initialization</h3>
01446 // The table array stores factors of n1 and n2, and trigonometric tables
01447 // that are used in calculation of the FFT.  table must be initialized by
01448 // calling the routine with isign = 0.  table does not have to be
01449 // reinitialized if the values of the problem sizes, n1 and n2, do not
01450 // change.
01451 //
01452 // <h3>Dimensions</h3>
01453 // In the preceding description, it is assumed that array subscripts were
01454 // zero-based, as is customary in FFT applications.  Thus, the input and
01455 // output arrays are declared:
01456 //
01457 // <srcblock>
01458 //      REAL    X(0:ldx-1, 0:n2-1)
01459 //      COMPLEX Y(0:ldy-1, 0:n2-1)
01460 // </srcblock>
01461 //
01462 // No change is made in the calling sequence, however, if you prefer to
01463 // use the more customary Fortran style with subscripts starting at 1.
01464 // The same values of ldx and ldy would be passed to the subroutine even
01465 // if the input and output arrays were dimensioned as follows:
01466 //
01467 // <srcblock>
01468 //      REAL    X(ldx, n2)
01469 //      COMPLEX Y(ldy, n2)
01470 // </srcblock>
01471 //
01472 // <example>
01473 // The following examples are for Origin series only.
01474 //
01475 // Example 1:  Initialize the TABLE array in preparation for doing a
01476 // two-dimensional FFT of size 128 by 256.  In this case, only the isign,
01477 // n1, n2, and table arguments are used; you can use dummy arguments or
01478 // zeros for other arguments.
01479 //
01480 // <srcblock>
01481 //       REAL TABLE ((15 + 128) + 2(15 + 256))
01482 //       CALL SCFFT2D (0, 128, 256, 0.0, DUMMY, 1, DUMMY, 1,
01483 //      &  TABLE, DUMMY, 0)
01484 // </srcblock>
01485 //
01486 // Example 2:  X is a real array of size (0:128, 0: 255), and Y is a
01487 // complex array of dimension (0:64, 0:255).  The first 128 elements of
01488 // each column of X contain data; for performance reasons, the extra
01489 // element forces the leading dimension to be an odd number.  Take the
01490 // two-dimensional FFT of X and store it in Y.  Initialize the TABLE
01491 // array, as in example 1.
01492 //
01493 // <srcblock>
01494 //       REAL    X(0:128, 0:255)
01495 //       COMPLEX Y(0:64, 0:255)
01496 //       REAL    TABLE ((15 + 128) + 2(15 + 256))
01497 //       REAL    WORK(128*256)
01498 //       ...
01499 //       CALL SCFFT2D(0, 128, 256, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
01500 //       CALL SCFFT2D(1, 128, 256, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
01501 // </srcblock>
01502 //
01503 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
01504 // and store it back in X.  The scale factor 1/(128*256) is used.  Assume
01505 // that the TABLE array is initialized already.
01506 //
01507 // <srcblock>
01508 //       CALL CSFFT2D(-1, 128, 256, 1.0/(128.0*256.0), Y, 65,
01509 //      &  X, 130, TABLE, WORK, 0)
01510 // </srcblock>
01511 //
01512 // Example 4:  Perform the same computation as in example 2, but assume
01513 // that the lower bound of each array is 1, rather than 0.  No change is
01514 // needed in the subroutine calls.
01515 //
01516 // <srcblock>
01517 //       REAL    X(129, 256)
01518 //       COMPLEX Y(65, 256)
01519 //       ...
01520 //       CALL SCFFT2D(0, 128, 256, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
01521 //       CALL SCFFT2D(1, 128, 256, 1.0, X, 129, Y, 65, TABLE, WORK, 0)
01522 // </srcblock>
01523 //
01524 // Example 5:  Perform the same computation as in example 4, but
01525 // equivalence the input and output arrays to save storage space.  In
01526 // this case, a row must be added to X, because it is equivalenced to a
01527 // complex array.  Assume that TABLE is already initialized.
01528 //
01529 // <srcblock>
01530 //       REAL    X(130, 256)
01531 //       COMPLEX Y(65, 256)
01532 //       EQUIVALENCE ( X(1, 1), Y(1, 1) )
01533 //       ...
01534 //       CALL SCFFT2D(1, 128, 256, 1.0, X, 130, Y, 65, TABLE, WORK, 0)
01535 // </srcblock>
01536 // </example>
01537 //
01538 // Input parameters:
01539 // <dl compact>
01540 // <dt><b>isign</b>
01541 // <dd>    Integer.
01542 //         Specifies whether to initialize the table array or to do the
01543 //         forward or inverse Fourier transform, as follows:
01544 //
01545 //         If isign = 0, the routine initializes the table array and
01546 //         returns.  In this case, the only arguments used or checked
01547 //         are isign, n, and table.
01548 //
01549 //         If isign = +1 or -1, the value of isign is the sign of the
01550 //         exponent used in the FFT formula.
01551 // <dt><b>n1</b>
01552 // <dd>    Integer.
01553 //         Transform size in the first dimension.  If n1 is not
01554 //         positive, <src>scfft2d</src> returns without calculating a
01555 //         transform.
01556 // <dt><b>n2</b>
01557 // <dd>    Integer.
01558 //         Transform size in the second dimension.  If n2 is not
01559 //         positive, <src>scfft2d</src> returns without calculating a
01560 //         transform.
01561 // <dt><b>scale</b>
01562 // <dd>    Scale factor.
01563 //         <src>scfft2d</src>: real. 
01564 //         <src>dzfft2d</src>: double precision.
01565 //         <src>csfft2d</src>: real.
01566 //         <src>zdfft2d</src>: double precision.
01567 //         Each element of the output array is multiplied by scale
01568 //         factor after taking the Fourier transform, as defined
01569 //         previously.
01570 // <dt><b>x</b>
01571 // <dd>    Array of dimension (0:ldx-1, 0:n2-1).
01572 //         <src>scfft2d</src>: real array.
01573 //         <src>dzfft2d</src>: double precision array.
01574 //         <src>csfft2d</src>: complex array.
01575 //         <src>zdfft2d</src>: double complex array.
01576 //
01577 //         Array of values to be transformed.
01578 // <dt><b>ldx</b>
01579 // <dd>    Integer.
01580 //         The number of rows in the x array, as it was declared in the
01581 //         calling program.  That is, the leading dimension of x.
01582 //         <src>scfft2d, dzfft2d</src>:  ldx >= MAX(n1, 1).
01583 //         <src>csfft2d, zdfft2d</src>:  ldx >= MAX(n1/2 + 1, 1).
01584 // <dt><b>ldy</b>
01585 // <dd>    Integer.
01586 //
01587 //         The number of rows in the y array, as it was declared in the
01588 //         calling program (the leading dimension of y).
01589 //
01590 //         <src>scfft2d, dzfft2d</src>:  ldy >= MAX(n1/2 + 1, 1).
01591 //         <src>csfft2d, zdfft2d</src>:  ldy >= MAX(n1 + 2, 1).
01592 //
01593 //         In the complex-to-real routine, two extra elements are in
01594 //         the first dimension (ldy >= n1 + 2, rather than just ldy >=
01595 //         n1).  These elements are needed for intermediate storage
01596 //         during the computation.  On exit, their value is undefined.
01597 // <dt><b>isys</b>
01598 // <dd>    Integer array of dimension (0:isys(0)).
01599 //         The first element of the array specifies how many more
01600 //         elements are in the array.  Use isys to specify certain
01601 //         processor-specific parameters or options.
01602 //
01603 //         If isys(0) = 0, the default values of such parameters are
01604 //         used.  In this case, you can specify the argument value as
01605 //         the scalar integer constant 0.
01606 //
01607 //         If isys(0) > 0, isys(0) gives the upper bound of the isys
01608 //         array.  Therefore, if il = isys(0), user-specified
01609 //         parameters are expected in isys(1) through isys(il).
01610 // <dt><b>isys</b>
01611 // <dd>    Algorithm used; value dependent on hardware system.  Currently, no
01612 //         special options are supported; therefore, you must always specify
01613 //         an isys argument as constant 0.
01614 // </dl>
01615 // Output parameters:
01616 // <dl compact>
01617 // <dt><b>y</b>
01618 // <dd>    <src>scfft2d</src>: complex array.
01619 //         <src>dzfft2d</src>: double complex array.
01620 //         <src>csfft2d</src>: real array.
01621 //         <src>zdfft2d</src>: double precision array.
01622 //
01623 //         Output array of transformed values.  The output array can be
01624 //         the same as the input array, in which case, the transform is
01625 //         done in place and the input array is overwritten with the
01626 //         transformed values.  In this case, it is necessary that the
01627 //         following equalities hold:
01628 //
01629 //         <src>scfft2d, dzfft2d</src>:  ldx = 2 * ldy.
01630 //         <src>csfft2d, zdfft2d</src>:  ldy = 2 * ldx.
01631 // <dt><b>table</b>
01632 // <dd>    Real array; dimension (15 + n1) + 2(15 + n2).
01633 //
01634 //         Table of factors and trigonometric functions.
01635 //
01636 //         If isign = 0, the routine initializes table (table is output
01637 //         only).
01638 //
01639 //         If isign = +1 or -1, the values in table are assumed to be
01640 //         initialized already by a prior call with isign = 0 (table is
01641 //         input only).
01642 // <dt><b>work</b>
01643 // <dd>    Real array; dimension (n1 * n2).
01644 //
01645 //         Work array.  This is a scratch array used for intermediate
01646 //         calculations.  Its address space must be different from that
01647 //         of the input and output arrays.
01648 // </dl>
01649 // <group>
01650 static void scfft2d(Int isign, Int n1, Int n2, Float scale, Float*
01651                     x, Int ldx, Complex* y, Int ldy, Float* table,
01652                     Float* work, Int isys); 
01653 static void dzfft2d(Int isign, Int n1, Int n2, Double scale, Double*
01654                     x, Int ldx, DComplex* y, Int ldy, Double* table,
01655                     Double* work, Int isys); 
01656 static void csfft2d(Int isign, Int n1, Int n2, Float scale, Complex*
01657                     x, Int ldx, Float* y, Int ldy, Float* table,
01658                     Float* work, Int isys); 
01659 static void zdfft2d(Int isign, Int n1, Int n2, Double scale, DComplex*
01660                     x, Int ldx, Double* y, Int ldy, Double* table,
01661                     Double* work, Int isys); 
01662 // </group>
01663 
01664 // These routines compute the three-dimensional complex FFT of the
01665 // complex matrix X, and store the results in the complex matrix Y.
01666 //
01667 // In FFT applications, it is customary to use zero-based subscripts; the
01668 // formulas are simpler that way.  So suppose the arrays are dimensioned
01669 // as follows:
01670 //
01671 // <srcblock>
01672 //      COMPLEX X(0:n1-1, 0:n2-1, 0:n3-1)
01673 //      COMPLEX Y(0:n1-1, 0:n2-1, 0:n3-1)
01674 // </srcblock>
01675 //
01676 // These routines compute the formula:
01677 //
01678 // <srcblock>
01679 //    Y(k1,k2,k3) =
01680 //        n1-1 n2-1 n3-1
01681 //    scale * Sum  Sum  Sum [X(j1,j2,j3)*w1**(j1*k1)*w2**(j2*k2)*w3**(j3*k3)]
01682 //        j1=0 j2=0 j3=0
01683 //
01684 //    for k1 = 0, ..., n1 - 1,
01685 //    k2 = 0, ..., n2 - 1,
01686 //    k3 = 0, ..., n3 - 1,
01687 //
01688 //    where:
01689 //    w1 = exp(isign*2*pi*i/n1),
01690 //    w2 = exp(isign*2*pi*i/n2),
01691 //    w3 = exp(isign*2*pi*i/n3),
01692 //    i  = + sqrt(-1)
01693 //    pi = 3.14159...
01694 //    isign = +1 or -1
01695 // </srcblock>
01696 //
01697 // Different authors use different conventions for which of the
01698 // transforms, isign = +1 or isign = -1, is the forward or inverse
01699 // transform, and what the scale factor should be in either case.  You
01700 // can make this routine compute any of the various possible definitions,
01701 // however, by choosing the appropriate values for isign and scale.
01702 //
01703 // The relevant fact from FFT theory is this:  If you take the FFT with
01704 // any particular values of isign and scale, the mathematical inverse
01705 // function is computed by taking the FFT with -isign and 1/(n1 * n2 * n3
01706 // * scale).  In particular, if you use isign = +1 and scale = 1.0 for
01707 // the forward FFT, you can compute the inverse FFT by using isign = -1
01708 // and scale = 1/(n1 . n2 . n3).
01709 //
01710 // <example>
01711 // The following examples are for Origin series only.
01712 //
01713 // Example 1:  Initialize the TABLE array in preparation for doing a
01714 // three-dimensional FFT of size 128 by 128 by 128.  In this case, only
01715 // the isign, n1, n2, n3, and table arguments are used; you can use dummy
01716 // arguments or zeros for other arguments.
01717 //
01718 // <srcblock>
01719 //       REAL TABLE ((30 + 256) + (30 + 256) + (30 + 256))
01720 //       CALL CCFFT3D (0, 128, 128, 128, 0.0, DUMMY, 1, 1, DUMMY, 1, 1,
01721 //      &  TABLE, DUMMY, 0)
01722 // </srcblock>
01723 //
01724 // Example 2:  X and Y are complex arrays of dimension (0:128, 0:128,
01725 // 0:128).  The first 128 elements of each dimension contain data; for
01726 // performance reasons, the extra element forces the leading dimensions
01727 // to be odd numbers.  Take the three-dimensional FFT of X and store it
01728 // in Y.  Initialize the TABLE array, as in example 1.
01729 //
01730 // <srcblock>
01731 //       COMPLEX X(0:128, 0:128, 0:128)
01732 //       COMPLEX Y(0:128, 0:128, 0:128)
01733 //       REAL TABLE ((30+256) + (30 + 256) + (30 + 256))
01734 //       REAL     WORK 2(128*128*128)
01735 //       ...
01736 //       CALL CCFFT3D(0, 128, 128, 128, 1.0, DUMMY, 1, 1,
01737 //     &    DUMMY, 1, 1, TABLE, WORK, 0)
01738 //       CALL CCFFT3D(1, 128, 128, 128, 1.0, X, 129, 129,
01739 //     &    Y, 129, 129, TABLE, WORK, 0)
01740 // </srcblock>
01741 //
01742 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
01743 // and store it back in X.  The scale factor 1.0/(128.0**3) is used.
01744 // Assume that the TABLE array is already initialized.
01745 //
01746 // <srcblock>
01747 //       CALL CCFFT3D(-1, 128, 128, 128, 1.0/(128.0**3), Y, 129, 129,
01748 //      &  X, 129, 129, TABLE, WORK, 0)
01749 // </srcblock>
01750 //
01751 // Example 4:  Perform the same computation as in example 2, but assume
01752 // that the lower bound of each array is 1, rather than 0.  The
01753 // subroutine calls do not change.
01754 //
01755 // <srcblock>
01756 //       COMPLEX X(129, 129, 129)
01757 //       COMPLEX Y(129, 129, 129)
01758 //       ...
01759 //       CALL CCFFT3D(0, 128, 128, 128, 1.0, DUMMY, 1, 1,
01760 //      &    DUMMY, 1, 1, TABLE, WORK, 0)
01761 //       CALL CCFFT3D(1, 128, 128, 128, 1.0, X, 129, 129,
01762 //      &    Y, 129, 129, TABLE, WORK, 0)
01763 // </srcblock>
01764 //
01765 // Example 5:  Perform the same computation as in example 4, but put the
01766 // output back in the array X to save storage space.  Assume that the
01767 // TABLE array is already initialized.
01768 //
01769 // <srcblock>
01770 //       COMPLEX X(129, 129, 129)
01771 //       ...
01772 //       CALL CCFFT3D(1, 128, 128, 128, 1.0, X, 129, 129,
01773 //      &    X, 129, 129, TABLE, WORK, 0)
01774 // </srcblock>
01775 // </example>
01776 //
01777 // Input parameters:
01778 // <dl compact>
01779 // <dt><b>isign</b>
01780 // <dd>    Integer.
01781 //         Specifies whether to initialize the table array or to do the
01782 //         forward or inverse Fourier transform, as follows:
01783 //
01784 //         If isign = 0, the routine initializes the table array and
01785 //         returns.  In this case, the only arguments used or checked
01786 //         are isign, n1, n2, n3, and table.
01787 //
01788 //         If isign = +1 or -1, the value of isign is the sign of the
01789 //         exponent used in the FFT formula.
01790 //
01791 // <dt><b>n1</b>
01792 // <dd>    Integer.
01793 //         Transform size in the first dimension.  If n1 is not
01794 //         positive, the routine returns without computing a transform.
01795 //
01796 // <dt><b>n2</b>
01797 // <dd>    Integer.
01798 //         Transform size in the second dimension.  If n2 is not
01799 //         positive, the routine returns without computing a transform.
01800 //
01801 // <dt><b>n3</b>
01802 // <dd>    Integer.
01803 //         Transform size in the third dimension.  If n3 is not
01804 //         positive, the routine returns without computing a transform.
01805 //
01806 // <dt><b>scale</b>
01807 // <dd>    Scale factor.
01808 //         <src>ccfft3d</src>: real.
01809 //         <src>zzfft3d</src>: double precision.
01810 //
01811 //         Each element of the output array is multiplied by scale
01812 //         after taking the Fourier transform, as defined previously.
01813 //
01814 // <dt><b>x</b>
01815 // <dd>    Array of dimension (0:ldx-1, 0:ldx2-1, 0:n3-1).
01816 //         <src>ccfft3d</src>: complex array.
01817 //         <src>zzfft3d</src>: double complex array.
01818 //
01819 //         Input array of values to be transformed.
01820 //
01821 // <dt><b>ldx</b>
01822 // <dd>    Integer.
01823 //         The first dimension of x, as it was declared in the calling
01824 //         program (the leading dimension of x).  ldx >= MAX(n1, 1).
01825 //
01826 // <dt><b>ldx2</b>
01827 // <dd>    Integer.
01828 //         The second dimension of x, as it was declared in the calling
01829 //         program.  ldx2 >= MAX(n2, 1).
01830 //
01831 // <dt><b>ldy</b>
01832 // <dd>    Integer.
01833 //         The first dimension of y, as it was declared in the calling
01834 //         program (the leading dimension of y).  ldy >= MAX(n1, 1).
01835 //
01836 // <dt><b>ldy2</b>
01837 // <dd>    Integer.
01838 //         The second dimension of y, as it was declared in the calling
01839 //         program.  ldy2 >= MAX(n2, 1).
01840 //
01841 // <dt><b>isys</b>
01842 // <dd>    Algorithm used; value dependent on hardware system.  Currently, no
01843 //         special options are supported; therefore, you must always specify
01844 //         an isys argument as constant 0.
01845 //
01846 //         isys = 0 or 1 depending on the amount of workspace the user
01847 //         can provide to the routine.
01848 // </dl>
01849 // Output parameters:
01850 // <dl compact>
01851 // <dt><b>y</b>
01852 // <dd>    Array of dimension (0:ldy-1, 0:ldy2-1, 0:n3-1).
01853 //         <src>ccfft3d</src>: complex array.
01854 //         <src>zzfft3d</src>: double complex array.
01855 //
01856 //         Output array of transformed values.  The output array may be
01857 //         the same as the input array, in which case, the transform is
01858 //         done in place; that is, the input array is overwritten with
01859 //         the transformed values.  In this case, it is necessary that
01860 //         ldx = ldy, and ldx2 = ldy2.
01861 //
01862 // <dt><b>table</b>
01863 // <dd>    Real array; dimension 30 + 2 * n1) + (30 + 2 * n2) + (30 + 2 * n3).
01864 //
01865 //         Table of factors and trigonometric functions.  If isign = 0,
01866 //         the routine initializes table (table is output only).  If
01867 //         isign = +1 or -1, the values in table are assumed to be
01868 //         initialized already by a prior call with isign = 0 (table is
01869 //         input only).
01870 //
01871 // <dt><b>work</b>
01872 // <dd>    Real array; dimension (n1 * n2 * n3).
01873 //
01874 //         Work array.  This is a scratch array used for intermediate
01875 //         calculations.  Its address space must be different from that
01876 //         of the input and output arrays.
01877 //
01878 // </dl>
01879 // <group>
01880 static void ccfft3d(Int isign, Int n1, Int n2, Int n3, Float scale,
01881                     Complex* x, Int ldx, Int ldx2, Complex* y, Int ldy,
01882                     Int ldy2, Float* table, Float* work, Int isys);
01883 static void zzfft3d(Int isign, Int n1, Int n2, Int n3, Double scale,
01884                     DComplex* x, Int ldx, Int ldx2, DComplex* y, Int
01885                     ldy, Int ldy2, Double* table, Double* work, Int
01886                     isys); 
01887 // </group>
01888 
01889 // These are C++ wrapper functions for the 3D real-to-complex and
01890 // complex-to-real transform routines in the SGI/Cray Scientific Library
01891 // (SCSL). The purpose of these definitions is to overload the functions so
01892 // that C++ users can access the functions in SCSL with identical function
01893 // names.
01894 //
01895 // <note role=warning> 
01896 // Currently, the SCSL is available only on SGI machines.
01897 // </note>
01898 //
01899 // <src>scfft3d/dzfft3d</src> computes the three-dimensional Fast Fourier
01900 // Transform (FFT) of the real matrix X, and it stores the results in the
01901 // complex matrix Y.  <src>csfft3d/zdfft3d</src> computes the corresponding
01902 // inverse transform.
01903 //
01904 // In FFT applications, it is customary to use zero-based subscripts; the
01905 // formulas are simpler that way.  First, the function of <src>SCFFT3D</src> is
01906 // described.  Suppose the arrays are dimensioned as follows:
01907 //
01908 // <srcblock>
01909 //      REAL    X(0:ldx-1, 0:ldx2-1, 0:n3-1)
01910 //      COMPLEX Y(0:ldy-1, 0:ldy2-1, 0:n3-1)
01911 // </srcblock>
01912 //
01913 // <src>scfft3d</src> computes the formula:
01914 //
01915 // <srcblock>
01916 //    Y(k1,k2,k3) =
01917 //        n1-1 n2-1 n3-1
01918 //    scale * Sum  Sum  Sum [X(j1,j2,j3)*w1**(j1*k1)*w2**(j2*k2)*w3**(j3*k3)]
01919 //        j1=0 j2=0 j3=0
01920 //
01921 //    for k1 = 0, ..., n1/2,
01922 //        k2 = 0, ..., n2 - 1,
01923 //        k3 = 0, ..., n3 - 1,
01924 //
01925 //    where:
01926 //        w1 = exp(isign*2*pi*i/n1),
01927 //        w2 = exp(isign*2*pi*i/n2),
01928 //        w3 = exp(isign*2*pi*i/n3),
01929 //        i = + sqrt(-1)
01930 //        pi = 3.14159...
01931 //        isign = +1 or -1
01932 // </srcblock>
01933 //
01934 // Different authors use different conventions for which of the
01935 // transforms, isign = +1 or isign = -1, is the forward or inverse
01936 // transform, and what the scale factor should be in either case.  You
01937 // can make these routines compute any of the various possible
01938 // definitions, however, by choosing the appropriate values for isign and
01939 // scale.
01940 //
01941 // The relevant fact from FFT theory is this:  If you take the FFT with
01942 // any particular values of isign and scale, the mathematical inverse
01943 // function is computed by taking the FFT with -isign and 1/(n1 * n2 * n3
01944 // * scale).  In particular, if you use isign = +1 and scale = 1.0 for
01945 // the forward FFT, you can compute the inverse FFT by isign = -1 and
01946 //
01947 // <srcblock>
01948 //      scale = 1.0/(n1*n2*n3).
01949 // </srcblock>
01950 //
01951 // <src>scfft3d</src> is very similar in function to <src>ccfft3d</src>, but
01952 // it takes the real-to-complex transform in the first dimension, followed by
01953 // the complex-to-complex transform in the second and third dimensions.
01954 //
01955 // <src>csfft3d</src> does the reverse.  It takes the complex-to-complex FFT
01956 // in the third and second dimensions, followed by the complex-to-real FFT in
01957 // the first dimension.
01958 //
01959 // See the <src>scfftm</src> man page for more information about
01960 // real-to-complex and complex-to-real FFTs.  The three dimensional analog of
01961 // the conjugate formula is as follows:
01962 //
01963 // <srcblock>
01964 //      Y         = conjg Y
01965 //       k ,k ,k            n1 - k , n2 - k , n3 - k
01966 //        1  2  3                 1        2        3
01967 //
01968 //      for  n1/2 <  k  <= n1 - 1
01969 //                    1
01970 //
01971 //           0 <= k  <= n2 - 1
01972 //                 2
01973 //
01974 //           0 <= k  <= n3 - 1
01975 //                 3
01976 // where the notation conjg(z) represents the complex conjugate of z.
01977 // </srcblock>
01978 //
01979 // Thus, you have to compute only (slightly more than) half out the
01980 // output values, namely:
01981 //
01982 // <srcblock>
01983 //      Y
01984 //       k ,k ,k
01985 //        1  2  3
01986 //
01987 //      for  0 <= k  <= n1/2
01988 //                 1
01989 //
01990 //           0 <= k  <= n2 - 1
01991 //                 2
01992 //
01993 //           0 <= k  <= n3 - 1
01994 // </srcblock>
01995 //
01996 // <h3>Algorithm</h3>
01997 // <src>scfft3d</src> uses a routine similar to <src>scfftm</src> to do
01998 // multiple FFTs first on all columns of the input matrix, then uses a routine
01999 // similar to <src>ccfftm</src> on all rows of the result, and then on all
02000 // planes of that result.  See <src>scfftm</src> and <src>ccfftm</src> for
02001 // more information about the algorithms used.
02002 //
02003 // </example>
02004 // The following examples are for Origin series only.
02005 //
02006 // Example 1:  Initialize the TABLE array in preparation for doing a
02007 // three-dimensional FFT of size 128 by 128 by 128.  In this case only
02008 // the isign, n1, n2, n3, and table arguments are used; you can use dummy
02009 // arguments or zeros for other arguments.
02010 //
02011 // <srcblock>
02012 //       REAL TABLE ((15 + 128) + 2(15+128) + 2( 15 + 128))
02013 //       CALL SCFFT3D (0, 128, 128, 128, 0.0, DUMMY, 1, 1, DUMMY, 1, 1,
02014 //      &  TABLE, DUMMY, 0)
02015 // </srcblock>
02016 //
02017 // Example 2:  X is a real array of size (0:128, 0:128, 0:128).  The
02018 // first 128 elements of each dimension contain data; for performance
02019 // reasons, the extra element forces the leading dimensions to be odd
02020 // numbers.  Y is a complex array of dimension (0:64, 0:128, 0:128).
02021 // Take the three-dimensional FFT of X and store it in Y.  Initialize the
02022 // TABLE array, as in example 1.
02023 //
02024 // <srcblock>
02025 //       REAL    X(0:128, 0:128, 0:128)
02026 //       COMPLEX Y(0:64,  0:128, 0:128)
02027 //       REAL TABLE ((15+128) + 2(15 + 128) + 2(15 + 128))
02028 //       REAL    WORK(128*128*128)
02029 //       ...
02030 //       CALL SCFFT3D(0, 128, 128, 128, 1.0, X, 129, 129,
02031 //      &    Y, 65, 129, TABLE, WORK, 0)
02032 //       CALL SCFFT3D(1, 128, 128, 128, 1.0, X, 129, 129,
02033 //      &    Y, 65, 129, TABLE, WORK, 0)
02034 // </srcblock>
02035 //
02036 // Example 3:  With X and Y as in example 2, take the inverse FFT of Y
02037 // and store it back in X.  The scale factor 1/(128**3) is used.  Assume
02038 // that the TABLE array is initialized already.
02039 //
02040 // <srcblock>
02041 //       CALL CSFFT3D(-1, 128, 128, 128, 1.0/128.0**3, Y, 65, 129,
02042 //      &    X, 130, 129, TABLE, WORK, 0)
02043 // </srcblock>
02044 //
02045 // Example 4:  Perform the same computation as in example 2, but assume
02046 // that the lower bound of each array is 1, rather than 0.  No change is
02047 // made in the subroutine calls.
02048 //
02049 // <srcblock>
02050 //       REAL    X(129, 129, 129)
02051 //       COMPLEX Y(65, 129, 129)
02052 //       REAL TABLE ((15+128) + 2(15 + 128) + 2(15 + 128))
02053 //       REAL    WORK(128*128*128)
02054 //       ...
02055 //       CALL SCFFT3D(0, 128, 128, 128, 1.0, X, 129, 129,
02056 //      &    Y, 65, 129, TABLE, WORK, 0)
02057 //       CALL SCFFT3D(1, 128, 128, 128, 1.0, X, 129, 129,
02058 //      &    X, 129, 129, TABLE, WORK, 0)
02059 // </srcblock>
02060 //
02061 // Example 5:  Perform the same computation as in example 4, but
02062 // equivalence the input and output arrays to save storage space.  Assume
02063 // that the TABLE array is initialized already.
02064 //
02065 // <srcblock>
02066 //       REAL    X(130, 129, 129)
02067 //       COMPLEX Y(65, 129, 129)
02068 //       EQUIVALENCE (X(1, 1, 1), Y(1, 1, 1))
02069 //       ...
02070 //       CALL SCFFT3D(1, 128, 128, 128, 1.0, X, 130, 129,
02071 //      &    Y, 65, 129, TABLE, WORK, 0)
02072 // </srcblock>
02073 // </example>
02074 //
02075 // Input parameters:
02076 // <dl compact>
02077 // <dt><b>isign</b>
02078 // <dd>    Integer.
02079 //         Specifies whether to initialize the table array or to do the
02080 //         forward or inverse Fourier transform, as follows:
02081 //
02082 //         If isign = 0, the routine initializes the table array and
02083 //         returns.  In this case, the only arguments used or checked
02084 //         are isign, n1, n2, n3, and table.
02085 //
02086 //         If isign = +1 or -1, the value of isign is the sign of the
02087 //         exponent used in the FFT formula.
02088 //
02089 // <dt><b>n1</b>
02090 // <dd>    Integer.
02091 //         Transform size in the first dimension.  If n1 is not
02092 //         positive, <src>scfft3d</src> returns without computing a transform.
02093 //
02094 // <dt><b>n2</b>
02095 // <dd>    Integer.
02096 //         Transform size in the second dimension.  If n2 is not
02097 //         positive, <src>scfft3d</src> returns without computing a transform.
02098 //
02099 // <dt><b>n3</b>
02100 // <dd>    Integer.
02101 //         Transform size in the third dimension.  If n3 is not
02102 //         positive, <src>scfft3d</src> returns without computing a transform.
02103 //
02104 // <dt><b>scale</b>
02105 // <dd>    Scale factor.
02106 //         <src>scfft3d</src>: real.
02107 //         <src>dzfft3d</src>: double precision.
02108 //         <src>csfft3d</src>: real.
02109 //         <src>zdfft3d</src>: double precision.
02110 //         Each element of the output array is multiplied by scale
02111 //         after taking the Fourier transform, as defined previously.
02112 //
02113 // <dt><b>x</b>
02114 // <dd>    Array of dimension (0:ldx-1, 0:ldx2-1, 0:n3-1).
02115 //         <src>scfft3d</src>: real array.
02116 //         <src>dzfft3d</src>: double precision array.
02117 //         <src>csfft3d</src>: complex array.
02118 //         <src>zdfft3d</src>: double complex array.
02119 //
02120 //         Array of values to be transformed.
02121 //
02122 // <dt><b>ldx</b>
02123 // <dd>    Integer.
02124 //         The first dimension of x, as it was declared in the calling
02125 //         program (the leading dimension of x).
02126 //
02127 //         <src>scfft3d, dzfft3d</src>:  ldx >= MAX(n1, 1).
02128 //         <src>csfft3d, zdfft3d</src>:  ldx >= MAX(n1/2 + 1, 1).
02129 //
02130 // <dt><b>ldx2</b>
02131 // <dd>    Integer.
02132 //         The second dimension of x, as it was declared in the calling
02133 //         program.  ldx2 >= MAX(n2, 1).
02134 //
02135 // <dt><b>ldy</b>
02136 // <dd>    Integer.
02137 //         The first dimension of y, as it was declared in the calling
02138 //         program; that is, the leading dimension of y.
02139 //
02140 //         <src>scfft3d, dzfft3d</src>:  ldy >= MAX(n1/2 + 1, 1).
02141 //         <src>csfft3d, zdfft3d</src>:  ldy >= MAX(n1 + 2, 1).
02142 //
02143 //         In the complex-to-real routine, two extra elements are in
02144 //         the first dimension (that is, ldy >= n1 + 2, rather than
02145 //         just ldy >= n1).  These elements are needed for intermediate
02146 //         storage during the computation.  On exit, their value is
02147 //         undefined.
02148 //
02149 // <dt><b>ldy2</b>
02150 // <dd>    Integer.
02151 //         The second dimension of y, as it was declared in the calling
02152 //         program.  ldy2 >= MAX(n2, 1).
02153 //
02154 // <dt><b>isys</b>
02155 // <dd>    Algorithm used; value dependent on hardware system.  Currently, no
02156 //         special options are supported; therefore, you must always specify
02157 //         an isys argument as constant 0.
02158 //
02159 //         isys = 0 or 1 depending on the amount of workspace the user
02160 //         can provide to the routine.
02161 // </dl>
02162 // Output parameters:
02163 // <dl compact>
02164 // <dt><b>y</b>
02165 // <dd>    Array of dimension (0:ldy-1, 0:ldy2-1, 0:n3-1).
02166 //         <src>scfft3d</src>: complex array.
02167 //         <src>dzfft3d</src>: double complex array.
02168 //         <src>csfft3d</src>: real array.
02169 //         <src>zdfft3d</src>: double precision array.
02170 //
02171 //         Output array of transformed values.  The output array can be
02172 //         the same as the input array, in which case, the transform is
02173 //         done in place; that is, the input array is overwritten with
02174 //         the transformed values.  In this case, it is necessary that
02175 //         the following equalities hold:
02176 //
02177 //         <src>scfft3d, dzfft3d</src>:  ldx = 2 * ldy, and ldx2 = ldy2.
02178 //         <src>csfft3d, zdfft3d</src>:  ldy = 2 * ldx, and ldx2 = ldy2.
02179 //
02180 // <dt><b>table</b>
02181 // <dd>    Real array; dimension (15 + n1) + 2(15 + n2) + 2(15 + n3).
02182 //
02183 //         Table of factors and trigonometric functions.
02184 //
02185 //         This array must be initialized by a call to <src>scfft3d</src> or
02186 //         <src>csfft3d</src> with isign = 0.
02187 //
02188 //         If isign = 0, table is initialized to contain trigonometric
02189 //         tables needed to compute a three-dimensional FFT of size n1
02190 //         by n2 by n3.  If isign = +1 or -1, the values in table are
02191 //         assumed to be initialized already by a prior call with isign
02192 //         = 0.
02193 //
02194 // <dt><b>work</b>
02195 // <dd>    Real array; dimension n1 * n2 * n3.
02196 //
02197 //         Work array.  This is a scratch array used for intermediate
02198 //         calculations.  Its address space must be different from that
02199 //         of the input and output arrays.
02200 //
02201 // </dl>
02202 // <group>
02203 static void scfft3d(Int isign, Int n1, Int n2, Int n3, Float scale,
02204                     Float* x, Int ldx, Int ldx2, Complex* y, Int ldy,
02205                     Int ldy2, Float* table, Float* work, Int isys); 
02206 static void dzfft3d(Int isign, Int n1, Int n2, Int n3, Double scale,
02207                     Double* x, Int ldx, Int ldx2, DComplex* y, Int
02208                     ldy, Int ldy2, Double* table, Double* work, Int
02209                     isys); 
02210 static void csfft3d(Int isign, Int n1, Int n2, Int n3, Float scale,
02211                     Complex* x, Int ldx, Int ldx2, Float* y, Int ldy,
02212                     Int ldy2, Float* table, Float* work, Int isys); 
02213 static void zdfft3d(Int isign, Int n1, Int n2, Int n3, Double scale,
02214                     DComplex* x, Int ldx, Int ldx2, Double* y, Int
02215                     ldy, Int ldy2, Double* table, Double* work, Int
02216                     isys); 
02217 // </group>
02218 };
02219 
02220 } //# NAMESPACE CASA - END
02221 
02222 #endif