FourierTransform.C 40.2 KB
Newer Older
Francois Gygi committed
1 2
////////////////////////////////////////////////////////////////////////////////
//
Francois Gygi committed
3 4 5 6
// Copyright (c) 2008 The Regents of the University of California
//
// This file is part of Qbox
//
Francois Gygi committed
7 8
// Qbox is distributed under the terms of the GNU General Public License
// as published by the Free Software Foundation, either version 2 of
Francois Gygi committed
9 10 11 12 13 14
// the License, or (at your option) any later version.
// See the file COPYING in the root directory of this distribution
// or <http://www.gnu.org/licenses/>.
//
////////////////////////////////////////////////////////////////////////////////
//
Francois Gygi committed
15 16 17
// FourierTransform.C
//
////////////////////////////////////////////////////////////////////////////////
18
// $Id: FourierTransform.C,v 1.21 2009-11-30 02:33:27 fgygi Exp $
19 20

// The following macros must be defined: USE_FFTW, USE_ESSL, USE_ESSL_2DFFT
Francois Gygi committed
21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36

#include "FourierTransform.h"
#include "Basis.h"
#include "Context.h"

#include <complex>
#include <algorithm>
#include <map>
#include <cassert>

#if USE_MPI
#include <mpi.h>
#else
typedef int MPI_Comm;
#endif

37
#if USE_FFTW
Francois Gygi committed
38
#include "fftw.h"
Francois Gygi committed
39 40 41
#ifdef ADD_
#define zdscal zdscal_
#endif
42 43

extern "C" void zdscal(int *n,double *alpha,std::complex<double> *x,int *incx);
44
#elif USE_ESSL
Francois Gygi committed
45
extern "C" {
46 47
  void dcft_(int *initflag, std::complex<double> *x, int *inc2x, int *inc3x,
             std::complex<double> *y, int *inc2y, int *inc3y,
Francois Gygi committed
48 49 50
             int *length, int *ntrans, int *isign,
             double *scale, double *aux1, int *naux1,
             double *aux2, int *naux2);
51 52
  void dcft2_(int *initflag, std::complex<double> *x, int *inc1x, int *inc2x,
             std::complex<double> *y, int *inc1y, int *inc2y,
Francois Gygi committed
53 54 55 56 57 58
             int *n1, int *n2, int *isign,
             double *scale, double *aux1, int *naux1,
             double *aux2, int *naux2);
#define USE_GATHER_SCATTER 1
}
#else
59
#define NO_FFT_LIB 1
60
void cfftm ( std::complex<double> *ain, std::complex<double> *aout,
61
  double scale, int ntrans, int length, int ainc, int ajmp, int idir );
Francois Gygi committed
62 63 64 65 66
#endif

#if USE_GATHER_SCATTER
extern "C" {
  // zgthr: x(i) = y(indx(i))
67
  void zgthr_(int* n, std::complex<double>* y,
68
              std::complex<double>* x, int*indx);
Francois Gygi committed
69
  // zsctr: y(indx(i)) = x(i)
70
  void zsctr_(int* n, std::complex<double>* x, int* indx,
71
              std::complex<double>* y);
Francois Gygi committed
72 73 74
}
#endif

75 76
using namespace std;

Francois Gygi committed
77 78 79
////////////////////////////////////////////////////////////////////////////////
FourierTransform::~FourierTransform()
{
80
#if USE_FFTW
Francois Gygi committed
81 82 83 84 85 86 87 88 89 90 91
  fftw_destroy_plan(fwplan0);
  fftw_destroy_plan(fwplan1);
  fftw_destroy_plan(fwplan2);
  fftw_destroy_plan(bwplan0);
  fftw_destroy_plan(bwplan1);
  fftw_destroy_plan(bwplan2);
#endif
}

////////////////////////////////////////////////////////////////////////////////
FourierTransform::FourierTransform (const Basis &basis,
92
  int np0, int np1, int np2) : ctxt_(basis.context()), basis_(basis),
Francois Gygi committed
93 94
  np0_(np0), np1_(np1), np2_(np2)
{
95
  assert(ctxt_.npcol() == 1);
Francois Gygi committed
96 97 98 99 100 101 102 103 104 105
  nprocs_ = ctxt_.size();
  myproc_ = ctxt_.myproc();

  //if ( ctxt_.onpe0() )
  //  cout << " FourierTransform: " << np0 << " " << np1 << " " << np2 << endl;

  np2_loc_.resize(nprocs_);
  np2_first_.resize(nprocs_);

  // Block-cyclic distribution for np2
106
  // Partition np2 into nprocs_ intervals and
Francois Gygi committed
107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134
  // store local sizes in np2_loc_[iproc]
  // Use same block distribution as in ScaLAPACK
  // Blocks 0,...,nprocs_-2 have size np2_block_size
  // Block nprocs_-1 may have a smaller size
  if ( np2_ % nprocs_ == 0 )
  {
    // all blocks have equal size
    const int np2_block_size = np2_ / nprocs_;
    for ( int iproc = 0; iproc < nprocs_; iproc++ )
      np2_loc_[iproc] = np2_block_size;
  }
  else
  {
    // first k-1 blocks have same size, k_th block is smaller, others zero
    const int np2_block_size = np2_ / nprocs_ + 1;
    const int k = np2_ / np2_block_size;
    for ( int iproc = 0; iproc < k; iproc++ )
      np2_loc_[iproc] = np2_block_size;
    np2_loc_[k] = np2_ - k * np2_block_size;
    for ( int iproc = k+1; iproc < nprocs_; iproc++ )
      np2_loc_[iproc] = 0;
  }

  np2_first_[0] = 0;
  for ( int iproc = 1; iproc < nprocs_; iproc++ )
  {
    np2_first_[iproc] = np2_first_[iproc-1] + np2_loc_[iproc-1];
  }
135

Francois Gygi committed
136 137 138 139 140 141 142 143 144 145 146 147 148
  // number of local z vectors
  if ( basis_.real() )
  {
    if ( myproc_ == 0 )
      // rod(0,0) is mapped to only one z vector
      nvec_ = 2 * basis_.nrod_loc() - 1;
    else
      nvec_ = 2 * basis_.nrod_loc();
  }
  else
  {
    nvec_ = basis_.nrod_loc();
  }
149

150 151 152 153 154 155
  // compute number of transforms along the x,y and z directions
  // ntrans0_ is the number of transforms along x in one of the two blocks
  // of vectors corresponding to positive and negative y indices
  ntrans0_ = max(abs(basis_.idxmax(1)),abs(basis_.idxmin(1)))+1;
  ntrans1_ = np0_;
  ntrans2_ = nvec_;
156

Francois Gygi committed
157 158
  // resize array zvec holding columns
  zvec_.resize(nvec_ * np2_);
159

Francois Gygi committed
160 161
  // Initialize FT library auxiliary arrays
  init_lib();
162

Francois Gygi committed
163 164
  // allocate send buffer
  sbuf.resize(nvec_ * np2_);
165

Francois Gygi committed
166 167 168 169 170
  // allocate receive buffer
  if ( basis_.real() )
    rbuf.resize((2 * basis_.nrods() - 1) * np2_loc_[myproc_]);
  else
    rbuf.resize(basis_.nrods() * np2_loc_[myproc_]);
171

Francois Gygi committed
172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205
  // compute send/receive counts and displacements in units of sizeof(double)

  scounts.resize(nprocs_);
  sdispl.resize(nprocs_);
  rcounts.resize(nprocs_);
  rdispl.resize(nprocs_);

  if ( basis_.real() )
  {
    for ( int iproc = 0; iproc < nprocs_; iproc++ )
    {
      scounts[iproc] = 2 * nvec_ * np2_loc_[iproc];
      int nvec_iproc = iproc == 0 ? 2*basis_.nrod_loc(iproc)-1 :
                                2 * basis_.nrod_loc(iproc);
      rcounts[iproc] = 2 * nvec_iproc * np2_loc_[myproc_];
    }
  }
  else
  {
    for ( int iproc = 0; iproc < nprocs_; iproc++ )
    {
      scounts[iproc] = 2 * nvec_ * np2_loc_[iproc];
      int nvec_iproc = basis_.nrod_loc(iproc);
      rcounts[iproc] = 2 * nvec_iproc * np2_loc_[myproc_];
    }
  }

  sdispl[0] = 0;
  rdispl[0] = 0;
  for ( int iproc = 1; iproc < nprocs_; iproc++ )
  {
    sdispl[iproc] = sdispl[iproc-1] + scounts[iproc-1];
    rdispl[iproc] = rdispl[iproc-1] + rcounts[iproc-1];
  }
206

Francois Gygi committed
207 208 209 210 211
  if ( basis_.real() )
  {
    // compute index arrays ifftp_ and ifftm_ for mapping vector->zvec
    ifftp_.resize(basis_.localsize());
    ifftm_.resize(basis_.localsize());
212

Francois Gygi committed
213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296
    if ( myproc_ == 0 )
    {
      // this process holds rod(0,0)
      // nvec_ == 2 * nrod_loc - 1

      // map rod(0,0)
      // the positive segment of rod(0,0) maps onto the first half of
      // the first column of zvec_, and the negative segment maps onto
      // the second half
      int ig = 0;
      ifftp_[0] = 0;
      ifftm_[0] = 0;
      ig++;
      for ( int l = 1; l < basis_.rod_size(0); l++ )
      {
        ifftp_[ig] = l;
        ifftm_[ig] = np2_ - l;
        ig++;
      }

      // map other rods(h,k) on pe 0, h!=0, k!=0
      // rod(h,k) maps onto column (2*irod-1)*np2_ of zvec_, irod=1,..,nrods-1
      // rod(-h,-k) maps onto column (2*irod)*np2_ of zvec_, irod=1,..,nrods-1
      for ( int irod = 1; irod < basis_.nrod_loc(); irod++ )
      {
        const int rodsize = basis_.rod_size(irod);
        for ( int i = 0; i < rodsize; i++ )
        {
          const int l = i + basis_.rod_lmin(irod);
          int izp =  l;
          int izm = -l;
          if ( izp < 0 ) izp += np2_;
          if ( izm < 0 ) izm += np2_;
          ifftp_[ig] = ( 2 * irod - 1 ) * np2_ + izp;
          ifftm_[ig] = ( 2 * irod ) * np2_ + izm;
          ig++;
        }
      }
      assert(ig == basis_.localsize());
    }
    else
    {
      // this process does not hold rod(0,0)
      // map rods(h,k)
      // rod(h,k)   maps onto column (2*irod)*np2_ of zvec_, irod=0,..,nrods-1
      // rod(-h,-k) maps onto column (2*irod+1)*np2_ of zvec_, irod=0,..,nrods-1
      int ig = 0;
      for ( int irod = 0; irod < basis_.nrod_loc(); irod++ )
      {
        const int rodsize = basis_.rod_size(irod);
        for ( int i = 0; i < rodsize; i++ )
        {
          const int l = i + basis_.rod_lmin(irod);
          int izp =  l;
          int izm = -l;
          if ( izp < 0 ) izp += np2_;
          if ( izm < 0 ) izm += np2_;
          ifftp_[ig] = ( 2 * irod ) * np2_ + izp;
          ifftm_[ig] = ( 2 * irod + 1 ) * np2_ + izm;
          ig++;
        }
      }
      assert(ig == basis_.localsize());
    }

    // compute ipack index array
    // used in packing zvec_ into sbuf
    // sbuf[ipack_[i]] = zvec_[i]
    ipack_.resize(nvec_*np2_);
    int idest = 0;
    for ( int iproc = 0; iproc < nprocs_; iproc++ )
    {
      int isource = np2_first_[iproc];
      int sz = np2_loc_[iproc];
      for ( int ivec = 0; ivec < nvec_; ivec++ )
      {
        for ( int i = 0; i < sz; i++ )
        {
          ipack_[isource+i] = idest + i;
        }
        idest += sz;
        isource += np2_;
      }
    }
297

Francois Gygi committed
298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314
    // compute array iunpack
    // used in unpacking rbuf into val
    // val[iunpack[i]] = rbuf[i]

    // rbuf contains 2*_nrods-1 segments of size np2_loc[myproc]
    // the position of vector ivec in local rbuf[_nrods*np2_loc_] is
    // obtained from rod_h[iproc][irod], rod_k[irod][iproc]
    // compute iunpack[i], i = 0, .. , _nrods * np2_loc_
    iunpack_.resize((2*basis_.nrods()-1)*np2_loc_[myproc_]);

    // map rod(0,0)
    for ( int l = 0; l < np2_loc_[myproc_]; l++ )
    {
      iunpack_[l] = l * np0_ * np1_;
    }
    int isource_p = np2_loc_[myproc_];
    int isource_m = 2 * np2_loc_[myproc_];
315

Francois Gygi committed
316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357
    // all rods of pe 0
    for ( int irod = 1; irod < basis_.nrod_loc(0); irod++ )
    {
      // map rod(h,k) and rod(-h,-k) columns of zvec_

      // map rod(h,k)
      // find position of rod(h,k) in the slab
      int hp = basis_.rod_h(0,irod);
      int kp = basis_.rod_k(0,irod);
      if ( hp < 0 ) hp += np0_;
      if ( kp < 0 ) kp += np1_;

      int hm = -hp;
      int km = -kp;
      if ( hm < 0 ) hm += np0_;
      if ( km < 0 ) km += np1_;

      for ( int l = 0; l < np2_loc_[myproc_]; l++ )
      {
        int idest_p = hp + np0_ * ( kp + np1_ * l );
        iunpack_[isource_p+l] = idest_p;

        int idest_m = hm + np0_ * ( km + np1_ * l );
        iunpack_[isource_m+l] = idest_m;
      }
      isource_p += 2 * np2_loc_[myproc_];
      isource_m += 2 * np2_loc_[myproc_];
    }

    // pe's above pe0
    for ( int iproc = 1; iproc < nprocs_; iproc++ )
    {
      for ( int irod = 0; irod < basis_.nrod_loc(iproc); irod++ )
      {
        // map rod(h,k) and rod(-h,-k) columns of zvec_

        // map rod(h,k)
        // find position of rod(h,k) in the slab
        int hp = basis_.rod_h(iproc,irod);
        int kp = basis_.rod_k(iproc,irod);
        if ( hp < 0 ) hp += np0_;
        if ( kp < 0 ) kp += np1_;
358

Francois Gygi committed
359 360 361 362
        int hm = -hp;
        int km = -kp;
        if ( hm < 0 ) hm += np0_;
        if ( km < 0 ) km += np1_;
363

Francois Gygi committed
364 365 366 367
        for ( int l = 0; l < np2_loc_[myproc_]; l++ )
        {
          int idest_p = hp + np0_ * ( kp + np1_ * l );
          iunpack_[isource_p+l] = idest_p;
368

Francois Gygi committed
369 370 371 372 373 374 375 376 377 378 379 380 381 382
          int idest_m = hm + np0_ * ( km + np1_ * l );
          iunpack_[isource_m+l] = idest_m;
        }
        isource_p += 2 * np2_loc_[myproc_];
        isource_m += 2 * np2_loc_[myproc_];
      }
    }
  }
  else
  {
    // basis is complex
    // compute index array ifftp_ for mapping vector->zvec
    // Note: ifftm_ is not used
    ifftp_.resize(basis_.localsize());
383

Francois Gygi committed
384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419
    // map rods(h,k)
    // rod(h,k)   maps onto column irod*np2_ of zvec_, irod=0,..,nrods-1
    int ig = 0;
    for ( int irod = 0; irod < basis_.nrod_loc(); irod++ )
    {
      const int rodsize = basis_.rod_size(irod);
      for ( int i = 0; i < rodsize; i++ )
      {
        const int l = i + basis_.rod_lmin(irod);
        int iz =  l;
        if ( iz < 0 ) iz += np2_;
        ifftp_[ig] = irod * np2_ + iz;
        ig++;
      }
    }
    assert(ig == basis_.localsize());

    // compute ipack index array
    // used in packing zvec_ into sbuf
    // sbuf[ipack_[i]] = zvec_[i]
    ipack_.resize(nvec_*np2_);
    int idest = 0;
    for ( int iproc = 0; iproc < nprocs_; iproc++ )
    {
      int isource = np2_first_[iproc];
      int sz = np2_loc_[iproc];
      for ( int ivec = 0; ivec < nvec_; ivec++ )
      {
        for ( int i = 0; i < sz; i++ )
        {
          ipack_[isource+i] = idest + i;
        }
        idest += sz;
        isource += np2_;
      }
    }
420

Francois Gygi committed
421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441
    // compute array iunpack
    // used in unpacking rbuf into val
    // val[iunpack[i]] = rbuf[i]

    // rbuf contains _nrods segments of size np2_loc[mype]
    // the position of vector ivec in local rbuf[_nrods*np2_loc_] is
    // obtained from rod_h[iproc][irod], rod_k[irod][iproc]
    // compute iunpack[i], i = 0, .. , _nrods * np2_loc_
    iunpack_.resize(basis_.nrods()*np2_loc_[myproc_]);

    int isource = 0;
    for ( int iproc = 0; iproc < nprocs_; iproc++ )
    {
      for ( int irod = 0; irod < basis_.nrod_loc(iproc); irod++ )
      {
        // map rod(h,k)
        // find position of rod(h,k) in the slab
        int h = basis_.rod_h(iproc,irod);
        int k = basis_.rod_k(iproc,irod);
        if ( h < 0 ) h += np0_;
        if ( k < 0 ) k += np1_;
442

Francois Gygi committed
443 444 445 446
        for ( int l = 0; l < np2_loc_[myproc_]; l++ )
        {
          int idest = h + np0_ * ( k + np1_ * l );
          iunpack_[isource+l] = idest;
447

Francois Gygi committed
448 449 450 451 452
        }
        isource += np2_loc_[myproc_];
      }
    }
  }
453

Francois Gygi committed
454 455 456 457 458
  // for ( int ig = 0; ig < basis_.localsize(); ig++ )
  // {
  //   assert(ifftp_[ig] >= 0 && ifftp_[ig] < zvec_.size());
  //   assert(ifftm_[ig] >= 0 && ifftm_[ig] < zvec_.size());
  // }
459

Francois Gygi committed
460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476
#if USE_GATHER_SCATTER
  // shift index array by one for fortran ZGTHR and ZSCTR calls
  for ( int i = 0; i < iunpack_.size(); i++ )
  {
    iunpack_[i]++;
  }
  for ( int i = 0; i < ipack_.size(); i++ )
  {
    ipack_[i]++;
  }
#endif
}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::backward(const complex<double>* c, complex<double>* f)
{
#if TIMING
477
  tm_b_map.start();
Francois Gygi committed
478 479 480
#endif
  vector_to_zvec(c);
#if TIMING
481
  tm_b_map.stop();
Francois Gygi committed
482 483 484 485 486 487 488 489
#endif
  bwd(f);
}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::forward(complex<double>* f, complex<double>* c)
{
  fwd(f);
490 491 492
#if TIMING
  tm_f_map.start();
#endif
Francois Gygi committed
493
  zvec_to_vector(c);
494 495 496
#if TIMING
  tm_f_map.stop();
#endif
Francois Gygi committed
497 498 499
}

////////////////////////////////////////////////////////////////////////////////
500
void FourierTransform::backward(const complex<double>* c1,
Francois Gygi committed
501 502 503
                               const complex<double>* c2,
                               complex<double>* f)
{
504 505 506
#if TIMING
  tm_b_map.start();
#endif
Francois Gygi committed
507
  doublevector_to_zvec(c1,c2);
508 509 510
#if TIMING
  tm_b_map.stop();
#endif
Francois Gygi committed
511 512 513 514 515 516 517 518
  bwd(f);
}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::forward(complex<double>* f,
  complex<double>* c1, complex<double>* c2)
{
  fwd(f);
519 520 521
#if TIMING
  tm_f_map.start();
#endif
Francois Gygi committed
522
  zvec_to_doublevector(c1,c2);
523 524 525
#if TIMING
  tm_f_map.stop();
#endif
Francois Gygi committed
526
}
527

Francois Gygi committed
528 529 530 531 532 533 534
////////////////////////////////////////////////////////////////////////////////
void FourierTransform::bwd(complex<double>* val)
{
  // transform zvec along z, transpose and transform along x,y, store
  // result in val
  // The columns of zvec[nvec_ * np2_] contain the full vectors
  // to be transformed
535 536
  //
  // If the basis is real: Column (h,k) is followed by column (-h,-k),
Francois Gygi committed
537 538 539
  // except for (0,0)

#if TIMING
540
  tm_b_fft.start();
Francois Gygi committed
541 542
#endif

543
#if USE_ESSL
Francois Gygi committed
544 545
  int inc1 = 1, inc2 = np2_, ntrans = nvec_, isign = -1, initflag = 0;
  double scale = 1.0;
546

Francois Gygi committed
547 548 549
  dcft_(&initflag,&zvec_[0],&inc1,&inc2,&zvec_[0],&inc1,&inc2,&np2_,&ntrans,
       &isign,&scale,&aux1zb[0],&naux1z,&aux2[0],&naux2);

550
#elif USE_FFTW
551
   /*
Francois Gygi committed
552 553 554 555 556 557 558 559 560 561 562 563
    * void fftw(fftw_plan plan, int howmany,
    *    FFTW_COMPLEX *in, int istride, int idist,
    *    FFTW_COMPLEX *out, int ostride, int odist);
    */
  int ntrans, inc1, inc2;

  ntrans = nvec_;
  inc1 = 1;
  inc2 = np2_;
  fftw(bwplan2,ntrans,(FFTW_COMPLEX*)&zvec_[0],inc1,inc2,
                      (FFTW_COMPLEX*)0,0,0);
#else
564 565 566 567 568 569 570 571 572
  // No library
  /* Transform along z */
  int ntrans = nvec_;
  int length = np2_;
  int ainc   = 1;
  int ajmp   = np2_;
  double scale = 1.0;
  int idir = -1;
  cfftm ( &zvec_[0], &zvec_[0], scale, ntrans, length, ainc, ajmp, idir );
Francois Gygi committed
573
#endif
574

Francois Gygi committed
575
#if TIMING
576 577
  tm_b_fft.stop();
  tm_b_pack.start();
Francois Gygi committed
578
#endif
579

Francois Gygi committed
580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603
  // scatter zvec_ to sbuf for transpose
#if USE_GATHER_SCATTER
  // zsctr: y(indx(i)) = x(i)
  // void zsctr_(int* n, complex<double>* x, int* indx, complex<double>* y);
  {
    complex<double>* y = &sbuf[0];
    complex<double>* x = &zvec_[0];
    int n = zvec_.size();
    zsctr_(&n,x,&ipack_[0],y);
  }
#else
  const int zvec_size = zvec_.size();
  double* const ps = (double*) &sbuf[0];
  const double* const pz = (double*) &zvec_[0];
  for ( int i = 0; i < zvec_size; i++ )
  {
    // sbuf[ipack_[i]] = zvec_[i];
    const int ip = ipack_[i];
    const double a = pz[2*i];
    const double b = pz[2*i+1];
    ps[2*ip]   = a;
    ps[2*ip+1] = b;
  }
#endif
604

Francois Gygi committed
605
  // segments of z-vectors are now in sbuf
606

Francois Gygi committed
607
#if TIMING
608 609
  tm_b_pack.stop();
  tm_b_mpi.start();
Francois Gygi committed
610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625
#endif

  // transpose
#if USE_MPI
  int status = MPI_Alltoallv((double*)&sbuf[0],&scounts[0],&sdispl[0],
      MPI_DOUBLE,(double*)&rbuf[0],&rcounts[0],&rdispl[0],MPI_DOUBLE,
      ctxt_.comm());
  if ( status != 0 )
  {
    cout << " FourierTransform: status = " << status << endl;
    ctxt_.abort(2);
  }
#else
  assert(sbuf.size()==rbuf.size());
  rbuf = sbuf;
#endif
626

Francois Gygi committed
627
#if TIMING
628 629
  tm_b_mpi.stop();
  tm_b_zero.start();
Francois Gygi committed
630 631 632 633 634 635 636 637 638 639 640 641 642
#endif

  // copy from rbuf to val
  // scatter index array iunpack
  {
    const int len = np012loc();
    double* const pv = (double*) &val[0];
    for ( int i = 0; i < len; i++ )
    {
      pv[2*i]   = 0.0;
      pv[2*i+1] = 0.0;
    }
  }
643

Francois Gygi committed
644
#if TIMING
645 646
  tm_b_zero.stop();
  tm_b_unpack.start();
Francois Gygi committed
647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672
#endif

#if USE_GATHER_SCATTER
  // zsctr(n,x,indx,y): y(indx(i)) = x(i)
  {
    complex<double>* y = &val[0];
    complex<double>* x = &rbuf[0];
    int n = rbuf.size();
    zsctr_(&n,x,&iunpack_[0],y);
  }
#else
  {
    const int rbuf_size = rbuf.size();
    const double* const pr = (double*) &rbuf[0];
    double* const pv = (double*) &val[0];
    for ( int i = 0; i < rbuf_size; i++ )
    {
      // val[iunpack_[i]] = rbuf[i];
      const int iu = iunpack_[i];
      const double a = pr[2*i];
      const double b = pr[2*i+1];
      pv[2*iu]   = a;
      pv[2*iu+1] = b;
    }
  }
#endif
673

Francois Gygi committed
674
#if TIMING
675 676
  tm_b_unpack.stop();
  tm_b_fft.start();
Francois Gygi committed
677 678 679 680 681
#endif

  for ( int k = 0; k < np2_loc_[myproc_]; k++ )
  {
    // transform along x for non-zero vectors only
682
    // transform along x for y in [0,ntrans0_] and y in [np1-ntrans0_, np1-1]
683 684
#if USE_ESSL
#if USE_ESSL_2DFFT
Francois Gygi committed
685 686 687 688

  // use 2D FFT for x and y transforms
  int inc1, inc2, istart, isign = -1, initflag = 0;
  double scale = 1.0;
689

Francois Gygi committed
690
  // xy transform
691
  istart = k * np0_ * np1_;
Francois Gygi committed
692 693 694 695 696 697 698
  inc1 = 1; inc2 = np0_;
  dcft2_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
        &np0_,&np1_,&isign,&scale,&aux1xyb[0],&naux1xy,&aux2[0],&naux2);

#else

  // use multiple 1-d FFTs for x and y transforms
699

Francois Gygi committed
700 701 702
  int inc1, inc2, ntrans, istart, length, isign = -1, initflag = 0;
  double scale = 1.0;
  // transform only non-zero vectors along x
703 704
  // First block: positive y indices: [0,ntrans0_]
  ntrans = ntrans0_;
Francois Gygi committed
705 706 707 708 709 710
  inc1 = 1;
  inc2 = np0_;
  istart = k * np0_ * np1_;
  length = np0_;
  dcft_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
       &length,&ntrans,&isign,&scale,&aux1xb[0],&naux1x,&aux2[0],&naux2);
711

712
  // Second block: negative y indices: [np1-ntrans0_,np1-1]
Francois Gygi committed
713 714
  inc1 = 1;
  inc2 = np0_;
715
  istart = np0_ * ( (np1_-ntrans) + k * np1_ );
Francois Gygi committed
716 717 718
  length = np0_;
  dcft_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
       &length,&ntrans,&isign,&scale,&aux1xb[0],&naux1x,&aux2[0],&naux2);
719

Francois Gygi committed
720 721 722 723 724 725 726 727
  // transform along y for all values of x
  ntrans = np0_;
  inc1 = np0_;
  inc2 = 1;
  istart = k * np0_ * np1_;
  length = np1_;
  dcft_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
       &length,&ntrans,&isign,&scale,&aux1yb[0],&naux1y,&aux2[0],&naux2);
728
#endif
729
#elif USE_FFTW
730
    int inc1, inc2, istart;
Francois Gygi committed
731

732 733
    int ntrans = ntrans0_;
    // Transform first block along x: positive y indices
Francois Gygi committed
734 735
    inc1 = 1;
    inc2 = np0_;
736
    istart = k * np0_ * np1_;
Francois Gygi committed
737 738
    fftw(bwplan0,ntrans,(FFTW_COMPLEX*)&val[istart],inc1,inc2,
                        (FFTW_COMPLEX*)0,0,0);
739
    // Transform second block along x: negative y indices
Francois Gygi committed
740 741
    inc1 = 1;
    inc2 = np0_;
742
    istart = np0_ * ( (np1_-ntrans) + k * np1_ );
Francois Gygi committed
743 744
    fftw(bwplan0,ntrans,(FFTW_COMPLEX*)&val[istart],inc1,inc2,
                        (FFTW_COMPLEX*)0,0,0);
745

Francois Gygi committed
746 747 748 749
    // transform along y for all values of x
    ntrans = np0_;
    inc1 = np0_;
    inc2 = 1;
750
    istart = k * np0_ * np1_;
Francois Gygi committed
751 752 753
    fftw(bwplan1,ntrans,(FFTW_COMPLEX*)&val[istart],inc1,inc2,
                        (FFTW_COMPLEX*)0,0,0);
#else
754 755
    // No library
    // transform along x for non-zero elements
756 757
    // Transform first block along x: positive y indices
    int ntrans = ntrans0_;
758 759 760 761 762 763 764
    int istart = k * np0_ * np1_;
    int length = np0_;
    int ainc   = 1;
    int ajmp   = np0_;
    double scale = 1.0;
    int idir = -1;
    cfftm (&val[istart],&val[istart],scale,ntrans,length,ainc,ajmp,idir );
765

766 767
    // Transform second block along x: negative y indices
    istart = np0_ * ( (np1_-ntrans) + k * np1_ );
768
    cfftm (&val[istart],&val[istart],scale,ntrans,length,ainc,ajmp,idir );
769

770 771 772 773 774 775 776 777
    // transform along y for all values of x
    ntrans = np0_;
    istart = k * np0_ * np1_;
    length = np1_;
    ainc = np0_;
    ajmp = 1;
    cfftm (&val[istart],&val[istart],scale,ntrans,length,ainc,ajmp,idir );

Francois Gygi committed
778
#endif
779
  } // for k
780

Francois Gygi committed
781
#if TIMING
782
  tm_b_fft.stop();
Francois Gygi committed
783 784 785 786 787 788
#endif
}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::fwd(complex<double>* val)
{
789 790 791
#if TIMING
  tm_f_fft.start();
#endif
Francois Gygi committed
792 793 794
  for ( int k = 0; k < np2_loc_[myproc_]; k++ )
  {
    // transform along x for non-zero vectors only
795
    // transform along x for y in [0,ntrans0_] and y in [np1-ntrans0_, np1-1]
796 797
#if USE_ESSL
#if USE_ESSL_2DFFT
Francois Gygi committed
798 799 800 801

    // use 2D FFT for x and y transforms
    int inc1, inc2, istart, isign = 1, initflag = 0;
    double scale = 1.0;
802

Francois Gygi committed
803 804 805 806 807 808 809 810 811
    // xy transform
    istart = k * np0_ * np1_;
    inc1 = 1; inc2 = np0_;
    dcft2_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
          &np0_,&np1_,&isign,&scale,&aux1xyf[0],&naux1xy,&aux2[0],&naux2);

#else

    // use multiple 1-d FFTs for x and y transforms
812

Francois Gygi committed
813 814 815 816 817 818 819 820 821 822
    int inc1, inc2, ntrans, istart, length, isign = 1, initflag = 0;
    double scale = 1.0;
    // transform along y for all values of x
    ntrans = np0_;
    inc1 = np0_;
    inc2 = 1;
    istart = k * np0_ * np1_;
    length = np1_;
    dcft_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
         &length,&ntrans,&isign,&scale,&aux1yf[0],&naux1y,&aux2[0],&naux2);
823

Francois Gygi committed
824
    // transform only non-zero vectors along x
825
    ntrans = ntrans0_;
Francois Gygi committed
826 827 828 829 830 831
    inc1 = 1;
    inc2 = np0_;
    istart = k * np0_ * np1_;
    length = np0_;
    dcft_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
         &length,&ntrans,&isign,&scale,&aux1xf[0],&naux1x,&aux2[0],&naux2);
832

Francois Gygi committed
833 834
    inc1 = 1;
    inc2 = np0_;
835
    istart = np0_ * ( (np1_-ntrans) + k * np1_ );
Francois Gygi committed
836 837 838
    length = np0_;
    dcft_(&initflag,&val[istart],&inc1,&inc2,&val[istart],&inc1,&inc2,
         &length,&ntrans,&isign,&scale,&aux1xf[0],&naux1x,&aux2[0],&naux2);
839
#endif
840
#elif USE_FFTW
841
    int inc1, inc2, istart;
Francois Gygi committed
842 843

    // transform along y for all values of x
844
    int ntrans = np0_;
Francois Gygi committed
845 846
    inc1 = np0_;
    inc2 = 1;
847
    istart = k * np0_ * np1_;
Francois Gygi committed
848 849
    fftw(fwplan1,ntrans,(FFTW_COMPLEX*)&val[istart],inc1,inc2,
                        (FFTW_COMPLEX*)0,0,0);
850

851 852
    ntrans = ntrans0_;
    // Transform first block along x: positive y indices
Francois Gygi committed
853 854
    inc1 = 1;
    inc2 = np0_;
855
    istart = k * np0_ * np1_;
Francois Gygi committed
856 857
    fftw(fwplan0,ntrans,(FFTW_COMPLEX*)&val[istart],inc1,inc2,
                        (FFTW_COMPLEX*)0,0,0);
858
    // Transform second block along x: negative y indices
Francois Gygi committed
859 860
    inc1 = 1;
    inc2 = np0_;
861
    istart = np0_ * ( (np1_-ntrans) + k * np1_ );
Francois Gygi committed
862 863
    fftw(fwplan0,ntrans,(FFTW_COMPLEX*)&val[istart],inc1,inc2,
                        (FFTW_COMPLEX*)0,0,0);
864

Francois Gygi committed
865
#else
866 867 868 869 870 871 872 873 874 875
    // No library
    // transform along y for all values of x
    int ntrans = np0_;
    int istart = k * np0_ * np1_;
    int length = np1_;
    int ainc = np0_;
    int ajmp = 1;
    double scale = 1.0;
    int idir = 1;
    cfftm (&val[istart],&val[istart],scale,ntrans,length,ainc,ajmp,idir );
876

877
    // transform along x for non-zero elements
878
    ntrans = ntrans0_;
879 880 881 882 883
    istart = k * np0_ * np1_;
    length = np0_;
    ainc   = 1;
    ajmp   = np0_;
    cfftm (&val[istart],&val[istart],scale,ntrans,length,ainc,ajmp,idir );
884

885
    istart = np0_ * ( (np1_-ntrans) + k * np1_ );
886
    cfftm (&val[istart],&val[istart],scale,ntrans,length,ainc,ajmp,idir );
Francois Gygi committed
887
#endif
888
  } // for k
889

Francois Gygi committed
890
  // gather val into rbuf
891 892 893 894
#if TIMING
  tm_f_fft.stop();
  tm_f_pack.start();
#endif
Francois Gygi committed
895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918

#if USE_GATHER_SCATTER
  // zgthr: x(i) = y(indx(i))
  // void zgthr_(int* n, complex<double>* y, complex<double>* x, int*indx);
  {
    complex<double>* y = &val[0];
    complex<double>* x = &rbuf[0];
    int n = rbuf.size();
    zgthr_(&n,y,x,&iunpack_[0]);
  }
#else
  const int rbuf_size = rbuf.size();
  double* const pr = (double*) &rbuf[0];
  const double* const pv = (double*) &val[0];
  for ( int i = 0; i < rbuf_size; i++ )
  {
    // rbuf[i] = val[iunpack_[i]];
    const int iu = iunpack_[i];
    const double a = pv[2*iu];
    const double b = pv[2*iu+1];
    pr[2*i]   = a;
    pr[2*i+1] = b;
  }
#endif
919

Francois Gygi committed
920
  // transpose
921 922 923 924
#if TIMING
  tm_f_pack.stop();
  tm_f_mpi.start();
#endif
Francois Gygi committed
925 926 927 928 929 930 931
#if USE_MPI
  int status = MPI_Alltoallv((double*)&rbuf[0],&rcounts[0],&rdispl[0],
      MPI_DOUBLE,(double*)&sbuf[0],&scounts[0],&sdispl[0],MPI_DOUBLE,
      ctxt_.comm());
  assert ( status == 0 );
#else
  assert(sbuf.size()==rbuf.size());
932
  sbuf = rbuf;
Francois Gygi committed
933
#endif
934

Francois Gygi committed
935 936
  // segments of z-vectors are now in sbuf
  // gather sbuf into zvec_
937 938 939 940
#if TIMING
  tm_f_mpi.stop();
  tm_f_unpack.start();
#endif
941

Francois Gygi committed
942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966
#if USE_GATHER_SCATTER
  // zgthr: x(i) = y(indx(i))
  // void zgthr_(int* n, complex<double>* y, complex<double>* x, int*indx);
  {
    complex<double>* y = &sbuf[0];
    complex<double>* x = &zvec_[0];
    int n = zvec_.size();
    zgthr_(&n,y,x,&ipack_[0]);
  }
#else
  const int zvec_size = zvec_.size();
  const double* const ps = (double*) &sbuf[0];
  double* const pz = (double*) &zvec_[0];
  for ( int i = 0; i < zvec_size; i++ )
  {
    // zvec_[i] = sbuf[ipack_[i]];
    const int ip = ipack_[i];
    const double a = ps[2*ip];
    const double b = ps[2*ip+1];
    pz[2*i]   = a;
    pz[2*i+1] = b;
  }
#endif

  // transform along z
967 968 969 970
#if TIMING
  tm_f_unpack.stop();
  tm_f_fft.start();
#endif
971

972
#if USE_ESSL
Francois Gygi committed
973 974
  int inc1 = 1, inc2 = np2_, ntrans = nvec_, isign = 1, initflag = 0;
  double scale = 1.0 / (np0_ * np1_ * np2_);
975

Francois Gygi committed
976 977
  dcft_(&initflag,&zvec_[0],&inc1,&inc2,&zvec_[0],&inc1,&inc2,&np2_,&ntrans,
        &isign,&scale,&aux1zf[0],&naux1z,&aux2[0],&naux2);
978

979
#elif USE_FFTW
Francois Gygi committed
980 981 982 983 984 985 986 987 988 989 990 991 992 993
 /*
  * void fftw(fftw_plan plan, int howmany,
  *    FFTW_COMPLEX *in, int istride, int idist,
  *    FFTW_COMPLEX *out, int ostride, int odist);
  */
  int ntrans, inc1, inc2;

  ntrans = nvec_;
  inc1 = 1;
  inc2 = np2_;
  fftw(fwplan2,ntrans,(FFTW_COMPLEX*)&zvec_[0],inc1,inc2,
                      (FFTW_COMPLEX*)0,0,0);
  int len = zvec_.size();
  double fac = 1.0 / ( np0_ * np1_ * np2_ );
Francois Gygi committed
994
  zdscal(&len,&fac,&zvec_[0],&inc1);
Francois Gygi committed
995
#else
996 997 998 999 1000 1001 1002 1003 1004
  // No library
  /* Transform along z */
  int ntrans = nvec_;
  int length = np2_;
  int ainc   = 1;
  int ajmp   = np2_;
  double scale = 1.0 / ( np0_ * np1_ * np2_ );
  int idir = 1;
  cfftm ( &zvec_[0], &zvec_[0], scale, ntrans, length, ainc, ajmp, idir );
Francois Gygi committed
1005
#endif
1006 1007 1008 1009

#if TIMING
  tm_f_fft.stop();
#endif
Francois Gygi committed
1010 1011 1012 1013 1014 1015 1016
}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::init_lib(void)
{
  // initialization of FFT libs

1017
#if USE_ESSL
1018
#if USE_ESSL_2DFFT
Francois Gygi committed
1019 1020 1021 1022 1023 1024 1025
  // use 2D FFT for x and y transforms and 1D FFT for z transforms
  naux1xy = 40000 + 2.28 * (np0_+np1_);
  aux1xyf.resize(naux1xy);
  aux1xyb.resize(naux1xy);
  int r = max(np0_,np1_);
  int s = min(64,min(np0_,np1_));
  naux2 = 20000 + (2*r+256)*(s+2.28);
1026

Francois Gygi committed
1027 1028 1029
  naux1z = 20000 + 2.28 * np2_;
  aux1zf.resize(naux1z);
  aux1zb.resize(naux1z);
1030

Francois Gygi committed
1031 1032 1033 1034
  int ntrans2 = nvec_;
  int naux2z = 20000 + 2.28 * np2_ + (256 + 2*np2_)*min(64,ntrans2);
  naux2 = max( naux2, naux2z );
  aux2.resize(naux2);
1035

Francois Gygi committed
1036
  double scale = 1.0;
1037

Francois Gygi committed
1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048
  // initialize xy transforms
  int initflag = 1, inc1, inc2, isign = -1;
  inc1 = 1; inc2 = np0_;
  dcft2_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np0_,&np1_,
         &isign,&scale,&aux1xyb[0],&naux1xy,&aux2[0],&naux2);
  isign = 1;
  dcft2_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np0_,&np1_,
         &isign,&scale,&aux1xyf[0],&naux1xy,&aux2[0],&naux2);

  // initialize z transforms
  int ntrans = nvec_;
1049
  inc1 = 1; inc2 = np2_;
Francois Gygi committed
1050 1051 1052 1053 1054 1055 1056 1057
  isign = -1;
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np2_,&ntrans,
       &isign,&scale,&aux1zb[0],&naux1z,&aux2[0],&naux2);
  isign = 1; scale = 1.0 / ( np0_ * np1_ * np2_ );
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np2_,&ntrans,
       &isign,&scale,&aux1zf[0],&naux1z,&aux2[0],&naux2);
#else

1058

1059 1060 1061
  naux1x = (int) (20000 + 2.28 * np0_);
  naux1y = (int) (20000 + 2.28 * np1_);
  naux1z = (int) (20000 + 2.28 * np2_);
Francois Gygi committed
1062 1063 1064 1065 1066 1067 1068
  aux1xf.resize(naux1x);
  aux1yf.resize(naux1y);
  aux1zf.resize(naux1z);
  aux1xb.resize(naux1x);
  aux1yb.resize(naux1y);
  aux1zb.resize(naux1z);

Francois Gygi committed
1069
  int naux2x = (int) (20000 + 2.28 * np0_ + (256 + 2*np0_)*min(64,ntrans0_));
Francois Gygi committed
1070
  naux2 = naux2x;
Francois Gygi committed
1071
  int naux2y = (int) (20000 + 2.28 * np1_ + (256 + 2*np1_)*min(64,ntrans1_));
Francois Gygi committed
1072
  naux2 = max( naux2, naux2y );
Francois Gygi committed
1073
  int naux2z = (int) (20000 + 2.28 * np2_ + (256 + 2*np2_)*min(64,ntrans2_));
Francois Gygi committed
1074 1075
  naux2 = max( naux2, naux2z );
  aux2.resize(naux2);
1076

Francois Gygi committed
1077 1078 1079 1080 1081
  // initialize x, y and z transforms

  int initflag = 1, inc1, inc2, ntrans, isign;
  double scale = 1.0;
  complex<double> *p = 0;
1082

Francois Gygi committed
1083
  // x transforms
1084
  inc1 = 1; inc2 = np0_; ntrans = ntrans0_;
Francois Gygi committed
1085 1086 1087 1088 1089 1090
  isign = -1;
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np0_,&ntrans,
        &isign,&scale,&aux1xb[0],&naux1x,&aux2[0],&naux2);
  isign = 1;
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np0_,&ntrans,
        &isign,&scale,&aux1xf[0],&naux1x,&aux2[0],&naux2);
1091

Francois Gygi committed
1092
  // y transforms
1093
  inc1 = np0_; inc2 = 1; ntrans = ntrans1_;
Francois Gygi committed
1094 1095 1096 1097 1098 1099
  isign = -1;
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np1_,&ntrans,
        &isign,&scale,&aux1yb[0],&naux1y,&aux2[0],&naux2);
  isign = 1;
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np1_,&ntrans,
        &isign,&scale,&aux1yf[0],&naux1y,&aux2[0],&naux2);
1100

Francois Gygi committed
1101
  // z transforms
1102
  inc1 = 1; inc2 = np2_; ntrans = ntrans2_;
Francois Gygi committed
1103 1104 1105 1106 1107 1108 1109
  isign = -1;
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np2_,&ntrans,
        &isign,&scale,&aux1zb[0],&naux1z,&aux2[0],&naux2);
  isign = 1; scale = 1.0 / ( np0_ * np1_ * np2_ );
  dcft_(&initflag,p,&inc1,&inc2,p,&inc1,&inc2,&np2_,&ntrans,
        &isign,&scale,&aux1zf[0],&naux1z,&aux2[0],&naux2);

1110 1111
#endif
#elif USE_FFTW
Francois Gygi committed
1112 1113

#if FFTWMEASURE
1114
  // FFTW_MEASURE
Francois Gygi committed
1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149
  fwplan0 = fftw_create_plan(np0_,FFTW_FORWARD,FFTW_MEASURE|FFTW_IN_PLACE);
  fwplan1 = fftw_create_plan(np1_,FFTW_FORWARD,FFTW_MEASURE|FFTW_IN_PLACE);
  fwplan2 = fftw_create_plan(np2_,FFTW_FORWARD,FFTW_MEASURE|FFTW_IN_PLACE);
  bwplan0 = fftw_create_plan(np0_,FFTW_BACKWARD,FFTW_MEASURE|FFTW_IN_PLACE);
  bwplan1 = fftw_create_plan(np1_,FFTW_BACKWARD,FFTW_MEASURE|FFTW_IN_PLACE);
  bwplan2 = fftw_create_plan(np2_,FFTW_BACKWARD,FFTW_MEASURE|FFTW_IN_PLACE);
#else
  // FFTW_ESTIMATE
  fwplan0 = fftw_create_plan(np0_,FFTW_FORWARD,FFTW_ESTIMATE|FFTW_IN_PLACE);
  fwplan1 = fftw_create_plan(np1_,FFTW_FORWARD,FFTW_ESTIMATE|FFTW_IN_PLACE);
  fwplan2 = fftw_create_plan(np2_,FFTW_FORWARD,FFTW_ESTIMATE|FFTW_IN_PLACE);
  bwplan0 = fftw_create_plan(np0_,FFTW_BACKWARD,FFTW_ESTIMATE|FFTW_IN_PLACE);
  bwplan1 = fftw_create_plan(np1_,FFTW_BACKWARD,FFTW_ESTIMATE|FFTW_IN_PLACE);
  bwplan2 = fftw_create_plan(np2_,FFTW_BACKWARD,FFTW_ESTIMATE|FFTW_IN_PLACE);
#endif

#else
  /* no library */
#endif

}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::vector_to_zvec(const complex<double> *c)
{
  const int ng = basis_.localsize();
  const int zvec_size = zvec_.size();
  double* const pz = (double*) &zvec_[0];
  for ( int i = 0; i < zvec_size; i++ )
  {
    pz[2*i]   = 0.0;
    pz[2*i+1] = 0.0;
  }
  const double* const pc = (double*) &c[0];
  if ( basis_.real() )
1150
  {
Francois Gygi committed
1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230
    for ( int ig = 0; ig < ng; ig++ )
    {
      // zvec_[ifftp_[ig]] = c[ig];
      // zvec_[ifftm_[ig]] = conj(c[ig]);
      const double a = pc[2*ig];
      const double b = pc[2*ig+1];
      const int ip = ifftp_[ig];
      const int im = ifftm_[ig];
      pz[2*ip] = a;
      pz[2*ip+1] = b;
      pz[2*im] = a;
      pz[2*im+1] = -b;
    }
  }
  else
    for ( int ig = 0; ig < ng; ig++ )
    {
      // zvec_[ifftp_[ig]] = c[ig];
      const double a = pc[2*ig];
      const double b = pc[2*ig+1];
      const int ip = ifftp_[ig];
      pz[2*ip] = a;
      pz[2*ip+1] = b;
    }
}
////////////////////////////////////////////////////////////////////////////////
void FourierTransform::zvec_to_vector(complex<double> *c)
{
  const int ng = basis_.localsize();
  const double* const pz = (double*) &zvec_[0];
  double* const pc = (double*) &c[0];
  for ( int ig = 0; ig < ng; ig++ )
  {
    // c[ig] = zvec_[ifftp_[ig]];
    const int ip = ifftp_[ig];
    const double pz0 = pz[2*ip];
    const double pz1 = pz[2*ip+1];
    pc[2*ig]   = pz0;
    pc[2*ig+1] = pz1;
  }
}

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::doublevector_to_zvec(const complex<double> *c1,
  const complex<double> *c2)
{
  // Mapping of two real functions onto zvec
  assert(basis_.real());
  const int zvec_size = zvec_.size();
  double* const pz = (double*) &zvec_[0];
  for ( int i = 0; i < zvec_size; i++ )
  {
    pz[2*i] = 0.0;
    pz[2*i+1] = 0.0;
  }
  const int ng = basis_.localsize();
  const double* const pc1 = (double*) &c1[0];
  const double* const pc2 = (double*) &c2[0];
  for ( int ig = 0; ig < ng; ig++ )
  {
    // const double a = c1[ig].real();
    // const double b = c1[ig].imag();
    // const double c = c2[ig].real();
    // const double d = c2[ig].imag();
    // zvec_[ip] = complex<double>(a-d, b+c);
    // zvec_[im] = complex<double>(a+d, c-b);
    const double a = pc1[2*ig];
    const double b = pc1[2*ig+1];
    const double c = pc2[2*ig];
    const double d = pc2[2*ig+1];
    const int ip = ifftp_[ig];
    const int im = ifftm_[ig];
    pz[2*ip]   = a - d;
    pz[2*ip+1] = b + c;
    pz[2*im]   = a + d;
    pz[2*im+1] = c - b;
  }
}

////////////////////////////////////////////////////////////////////////////////
1231
void FourierTransform::zvec_to_doublevector(complex<double> *c1,
Francois Gygi committed
1232 1233
  complex<double> *c2 )
{
1234
  // Mapping of zvec onto two real functions
Francois Gygi committed
1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261
  assert(basis_.real());
  const int ng = basis_.localsize();
  const double* const pz = (double*) &zvec_[0];
  double* const pc1 = (double*) &c1[0];
  double* const pc2 = (double*) &c2[0];
  for ( int ig = 0; ig < ng; ig++ )
  {
    // const double a = 0.5*zvec_[ip].real();
    // const double b = 0.5*zvec_[ip].imag();
    // const double c = 0.5*zvec_[im].real();
    // const double d = 0.5*zvec_[im].imag();
    // c1[ig] = complex<double>(a+c, b-d);
    // c2[ig] = complex<double>(b+d, c-a);
    const int ip = ifftp_[ig];
    const int im = ifftm_[ig];
    const double a = pz[2*ip];
    const double b = pz[2*ip+1];
    const double c = pz[2*im];
    const double d = pz[2*im+1];
    pc1[2*ig]   = 0.5 * ( a + c );
    pc1[2*ig+1] = 0.5 * ( b - d );
    pc2[2*ig]   = 0.5 * ( b + d );
    pc2[2*ig+1] = 0.5 * ( c - a );
  }
}


1262
#if NO_FFT_LIB
Francois Gygi committed
1263 1264 1265 1266

////////////////////////////////////////////////////////////////////////////////
//
//     /* no library: use cfftm function */
1267 1268
//
//
Francois Gygi committed
1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286
//     /* Transform along x */
//     int i2;
//     int ntrans = np1*np2;
//     int length = np0;
//     int ainc   = 1;
//     int ajmp   = np0i;
//     int idir;
//     double scale;
//     if ( dir == R_TO_K )
//     {
//       idir = -1;
//       scale = 1.0 / ((double) np0*np1*np2);
//     }
//     else
//     {
//       idir = 1;
//       scale = 1.0;
//     }
1287
//
Francois Gygi committed
1288
//     cfftm ( &val[0], &val[0], scale, ntrans, length, ainc, ajmp, idir );
1289
//
Francois Gygi committed
1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300
//     /* Transform along y */
//     for ( i2 = 0; i2 < np2; i2++ )
//     {
//       int ist = i2 * np0i * np1;
//       ntrans = np0;
//       length = np1;
//       ainc   = np0i;
//       ajmp   = 1;
//       scale = 1.0;
//       cfftm ( &val[ist], &val[ist], scale, ntrans, length, ainc, ajmp, idir );
//     }
1301
//
Francois Gygi committed
1302 1303 1304 1305 1306 1307 1308
//     /* Transform along z */
//     ntrans = np0i*np1;
//     length = np2;
//     ainc   = np0i*np1;
//     ajmp   = 1;
//     scale = 1.0;
//     cfftm ( &val[0], &val[0], scale, ntrans, length, ainc, ajmp, idir );
1309
//
Francois Gygi committed
1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361
////////////////////////////////////////////////////////////////////////////////

/* define multiple FFT function here */

void cfftm ( complex<double> *ain, complex<double> *aout, double scale,
  int ntrans, int length,
  int ainc, int ajmp, int idir )
/*
 *  cfftm: multiple one-dimensional complex FFT
 *
 *  ain     complex array (input)
 *  aout    complex array (output)
 *  scale   global scaling factor
 *  ntrans  number of transforms
 *  length  length of each transform (in complex numbers)
 *  ainc    distance between elements within a transform (in complex numbers)
 *  ajmp    distance between first elements of transforms (in complex numbers)
 *  idir    direction of transform
 */

{
  void cfft ( int idir, complex<double> *z1, complex<double> *z2, int n,
    int *inzee );
  vector<complex<double> > tmpa(length), tmpb(length);
  for ( int it = 0; it < ntrans; it++ )
  {
    int ibase = it * ajmp;
    for ( int i = 0; i < length; i++ )
    {
      tmpa[i] = ain[ibase+i*ainc];
    }
    int inzee = 1;
    cfft ( idir, &tmpa[0], &tmpb[0], length, &inzee );
    if ( inzee == 1 )
      for ( int i = 0; i < length; i++ )
      {
        aout[ibase+i*ainc] = tmpa[i];
      }
    else
      for ( int i = 0; i < length; i++ )
      {
        aout[ibase+i*ainc] = tmpb[i];
      }
    for ( int i = 0; i < length; i++ )
    {
      aout[ibase+i*ainc] *= scale;
    }
  }
}

/*******************************************************************************
 *
1362 1363
 *  Complex FFT
 *  C version
Francois Gygi committed
1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379
 *
 *  From: C++ Language System Release 3.0 Library Manual
 *  Transcription from 'FFT as Nested Multiplication, with a twist'
 *  C. de Boor, SIAM Sci. Stat. Comput., Vol 1, No 1, March 1980
 *
 *  Adapted to C by F.Gygi, 17 Feb 1993, 9 Dec 1993
 *
 ******************************************************************************/

#include <math.h>

#define NEXTMX 12

void cfft ( int idir, complex<double> *z1, complex<double> *z2, int n,
  int *inzee )
{
1380
  // Compute the discrete Fourier transform of z1 (or z2) in
Francois Gygi committed
1381 1382 1383 1384
  // the Cooley-Tukey way, but with a twist.
  // z1[before], z2[before]
  // *inzee == 1 means input in z1; *inzee == 2 means input in z2

1385
  void fftstp ( int idir, complex<double> *zin, int after,
Francois Gygi committed
1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401
                int now, int before, complex<double> *zout );
  static int prime[NEXTMX] = { 2, 3, 5, 7, 11, 13, 17, 19, 23, 29, 31, 37 };
  int before = n;
  int after = 1;
  int next = 0;
  int now;


  do
  {
    int np = prime[next];
    if ( (before/np)*np < before )
    {
      if ( ++next < NEXTMX ) continue;
      now = before;
      before = 1;
1402
    }
Francois Gygi committed
1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417
    else
    {
      now = np;
      before /= np;
    }
    if ( *inzee == 1 )
      fftstp ( idir, z1, after, now, before, z2 );
    else
      fftstp ( idir, z2, after, now, before, z1 );
    *inzee = 3 - *inzee;
    after *= now;
  } while ( before > 1 );

}

1418
void fftstp ( int idir, complex<double> *zin, int after,
Francois Gygi committed
1419 1420 1421 1422 1423 1424
              int now, int before, complex<double> *zout )
{

  static const double twopi = 2 * 3.141592653589793;
  double angle;
  complex<double> omega;
1425
  complex<double> arg,value;
Francois Gygi committed
1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451
  int ia,ib,j,k,l,in;

  angle = twopi/(now*after);
  omega =  complex<double>(cos ( angle ),-idir * sin ( angle ));
  arg = 1.0;
  for ( j = 0; j < now; j++ )
  {
    for ( ia = 0; ia < after; ia++ )
    {
      for ( ib = 0; ib < before; ib++ )
      {
        /* value = zin(ia,ib,now) */
        k = ia + ib*after + (now-1)*before*after;
        value = zin[k];
        for ( in = now-2; in >= 0; in-- )
        {
          /* value = value*arg + zin(ia,ib,in) */
          /* zin(ia,ib,in) = zin[ia + ib*after + in*before*after]; */
          l = ia + ib*after + in*before*after;
          value = value * arg + zin[l];
        }
        /* zout(ia,j,ib) = value */
        /* zout[ia + j*after + ib*now*after] = value; */
        l = ia + j*after + ib*now*after;
        zout[l] = value;
      }
1452

Francois Gygi committed
1453 1454 1455 1456 1457 1458 1459
      /* arg *= omega; */
      arg *= omega;

    }
  }
}
#endif
1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479

////////////////////////////////////////////////////////////////////////////////
void FourierTransform::reset_timers(void)
{
#if TIMING
  tm_f_map.reset();
  tm_f_fft.reset();
  tm_f_pack.reset();
  tm_f_mpi.reset();
  tm_f_zero.reset();
  tm_f_unpack.reset();

  tm_b_map.reset();
  tm_b_fft.reset();
  tm_b_pack.reset();
  tm_b_mpi.reset();
  tm_b_zero.reset();
  tm_b_unpack.reset();
#endif
}