Bisection.C 22.1 KB
Newer Older
Francois Gygi committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
////////////////////////////////////////////////////////////////////////////////
//
// Copyright (c) 2011 The Regents of the University of California
//
// This file is part of Qbox
//
// Qbox is distributed under the terms of the GNU General Public License
// as published by the Free Software Foundation, either version 2 of
// the License, or (at your option) any later version.
// See the file COPYING in the root directory of this distribution
// or <http://www.gnu.org/licenses/>.
//
////////////////////////////////////////////////////////////////////////////////
//
// Bisection.C
//
////////////////////////////////////////////////////////////////////////////////
#include "Bisection.h"
#include <bitset>
#include <algorithm>
21 22 23 24
#include "jade.h"
#include "FourierTransform.h"

using namespace std;
Francois Gygi committed
25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76

////////////////////////////////////////////////////////////////////////////////
int walsh(int l, int n, int i)
{
  // walsh function at level l at position i in a grid of n points
  assert(i>=0);
  assert(i<n);
  if ( l == 0 )
  {
    assert(n%2==0);
    if ( i >= n/2 ) return 0;
    else return 1;
  }
  else if ( l == 1 )
  {
    assert(n%4==0);
    if ( i >= n/4 && i < (3*n)/4 ) return 0;
    else return 1;
  }
  else if ( l == 2 )
  {
    assert(n%8==0);
    if ( (i >= n/8 && i < 3*n/8) || (i >= 5*n/8 && i < 7*n/8) ) return 0;
    else return 1;
  }
  else if ( l == 3 )
  {
    assert(n%16==0);
    if ( (i >= n/16 && i < 3*n/16) ||
         (i >= 5*n/16 && i < 7*n/16) ||
         (i >= 9*n/16 && i < 11*n/16) ||
         (i >= 13*n/16 && i < 15*n/16) ) return 0;
    else return 1;
  }
  else if ( l == 4 )
  {
    assert(n%32==0);
    if ( (i >= n/32 && i < 3*n/32) ||
         (i >= 5*n/32 && i < 7*n/32) ||
         (i >= 9*n/32 && i < 11*n/32) ||
         (i >= 13*n/32 && i < 15*n/32) ||
         (i >= 17*n/32 && i < 19*n/32) ||
         (i >= 21*n/32 && i < 23*n/32) ||
         (i >= 25*n/32 && i < 27*n/32) ||
         (i >= 29*n/32 && i < 31*n/32) ) return 0;
    else return 1;
  }
  else
    assert(false);
}

////////////////////////////////////////////////////////////////////////////////
77
Bisection::Bisection(const SlaterDet& sd, int nlevels[3]) : ctxt_(sd.context())
Francois Gygi committed
78 79
{
  // localization vectors are long int
80
  assert(sizeof(long int) >= 4);
Francois Gygi committed
81

82 83
  nst_ = sd.nst();
  nstloc_ = sd.nstloc();
Francois Gygi committed
84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99

  // number of bisection levels in each direction
  nlevels_[0]=nlevels[0];
  nlevels_[1]=nlevels[1];
  nlevels_[2]=nlevels[2];

  // number of subdivisions required in each direction
  // ndiv = 2^nlevel
  ndiv_[0] = 1 << nlevels[0];
  ndiv_[1] = 1 << nlevels[1];
  ndiv_[2] = 1 << nlevels[2];

  // largest number of levels
  nlevelsmax_=max(nlevels[0],max(nlevels[1],nlevels[2]));

  // real-space grid size for wave functions
100 101 102
  np_[0] = sd.basis().np(0);
  np_[1] = sd.basis().np(1);
  np_[2] = sd.basis().np(2);
Francois Gygi committed
103 104 105 106 107 108 109 110 111
  // adapt the grid dimensions to the levels of bisection
  for ( int idim = 0; idim < 3; idim++ )
  {
    for ( int j=1; j<=nlevels[idim]; j++ )
    {
      int base = 1 << j;
      if ( np_[idim] % base != 0 ) np_[idim] += base/2;
    }
  }
112 113 114
  while (!sd.basis().factorizable(np_[0])) np_[0] += (1<<nlevels[0]);
  while (!sd.basis().factorizable(np_[1])) np_[1] += (1<<nlevels[1]);
  while (!sd.basis().factorizable(np_[2])) np_[2] += (1<<nlevels[2]);
Francois Gygi committed
115 116

  // number of grid points of augmented grid for normalization
117
  ft_ = new FourierTransform(sd.basis(),np_[0],np_[1],np_[2]);
Francois Gygi committed
118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137
  np01_ = np_[0]*np_[1];
  np2loc_ = ft_->np2_loc();
  np012loc_ = ft_->np012loc();

  // xy projector index: xy_proj[i+j*np0] is the xy projector
  // associated with a point located at i + j*np0 + k*np0*np1
  xy_proj_.resize(np01_);

  for ( int i = 0; i < np_[0]; i++ )
    for ( int j = 0; j < np_[1]; j++ )
    {
      int i_slice_x= ( i * ndiv_[0] ) / np_[0];
      int i_slice_y= ( j * ndiv_[1] ) / np_[1];
      int xy_proj_index = i_slice_x + ndiv_[0] * i_slice_y;
      xy_proj_[i+np_[0]*j] = xy_proj_index;
    }

  // nmat_: number of A matrices
  nmat_ = nlevels[0] + nlevels[1] + nlevels[2];

138 139 140 141 142
  // each projector uses two bits of the localization vector
  // check that sizeof(long int) is sufficient
  // number of bits in a long int is 8*sizeof(long int)
  assert(2*nmat_ <= 8*sizeof(long int));

Francois Gygi committed
143 144
  // allocate A matrices
  amat_.resize(nmat_);
145
  const ComplexMatrix &c = sd.c();
Francois Gygi committed
146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161
  for ( int i = 0; i < nmat_; i++ )
    amat_[i] = new DoubleMatrix(c.context(),c.n(),c.n(),c.nb(),c.nb());
  // allocate rotation matrix
  u_ = new DoubleMatrix(c.context(),c.n(),c.n(),c.nb(),c.nb());

  // matrices of real space wave functions in subdomains
  rmat_.resize( ndiv_[0]*ndiv_[1] );
  {
    // xyproj_rsize: real-space size of xy projectors
    vector<int> xyproj_rsize( ndiv_[0]*ndiv_[1] , 0 );
    for ( int ixy = 0; ixy < np01_; ixy++ )
      xyproj_rsize[xy_proj_[ixy]] += np2loc_;

    // max_xyproj_rsize: maximum real-space size of xy projectors
    vector<int> max_xyproj_rsize( ndiv_[0]*ndiv_[1] , 0 );
    MPI_Allreduce((void*)&xyproj_rsize[0] , (void*)&max_xyproj_rsize[0],
162
      (int)xyproj_rsize.size(), MPI_INT ,MPI_MAX , ctxt_.comm());
Francois Gygi committed
163 164 165 166

    // allocate matrices rmat_[i]
    for ( int i = 0; i < rmat_.size(); i++ )
    {
167 168 169
      int n = c.n();
      int nb = c.nb();
      int m = max_xyproj_rsize[i] * c.context().nprow();
Francois Gygi committed
170
      int mb = max_xyproj_rsize[i];
171
      rmat_[i] = new DoubleMatrix(c.context(),m,n,mb,nb);
Francois Gygi committed
172 173 174 175
      rmat_[i]->clear();
    }
  }

176
  localization_.resize(nst_);
Francois Gygi committed
177 178 179 180 181 182 183 184 185 186 187 188
}

////////////////////////////////////////////////////////////////////////////////
Bisection::~Bisection(void)
{
  delete ft_;
  for ( int i = 0; i < nmat_; i++ ) delete amat_[i];
  delete u_;
  for ( int i = 0; i < rmat_.size(); i++ ) delete rmat_[i];
}

////////////////////////////////////////////////////////////////////////////////
189
void Bisection::compute_transform(const SlaterDet& sd, int maxsweep, double tol)
Francois Gygi committed
190
{
191 192
  // compute the transformation with tolerance tol

Francois Gygi committed
193 194
  map<string,Timer> tmap;

195 196 197
  // check that basis is real
  // jade is not implemented for complex matrices
  assert( sd.basis().real() );
Francois Gygi committed
198

199 200 201
  const ComplexMatrix& c = sd.c();
  const vector<double>& occ = sd.occ();
  assert(occ.size() == nst_);
Francois Gygi committed
202 203

#ifdef TIMING
204
  tmap["wf FFT"].start();
Francois Gygi committed
205
#endif
206 207 208 209 210 211 212 213 214
  // Fourier transform states and save real-space values in
  // rmat_ matrices
  vector<complex<double> > wftmp(np012loc_,0.0);
  for ( int n = 0; n < nstloc_; n++ )
  {
    ft_->backward(c.cvalptr(c.mloc()*n),&wftmp[0]);
    // pointers to rmat
    vector<double *> p_rmat( ndiv_[0]*ndiv_[1] );
    for (int iproj=0; iproj<ndiv_[0]*ndiv_[1]; iproj++)
Francois Gygi committed
215
    {
216 217 218 219 220 221 222
      int index = n * rmat_[iproj]->mloc();
      p_rmat[iproj] = rmat_[iproj]->valptr(index);
    }
    // copy wf to rmat arrays
    for ( int iz = 0; iz < np2loc_; iz++ )
    {
      for ( int ixy = 0; ixy < np01_; ixy++ )
Francois Gygi committed
223
      {
224 225 226
        int xy_proj = xy_proj_[ixy];
        (*p_rmat[xy_proj]) = wftmp[ixy + iz*np01_].real();
        p_rmat[xy_proj]++;
Francois Gygi committed
227 228
      }
    }
229
  }
Francois Gygi committed
230
#ifdef TIMING
231
  tmap["wf FFT"].stop();
Francois Gygi committed
232 233
#endif

234
  // compute the x/y A matrices
Francois Gygi committed
235
#ifdef TIMING
236
  tmap["xy products"].start();
Francois Gygi committed
237
#endif
238 239 240
  // clear a matrices
  for ( int i = 0; i < nmat_; i++ )
    amat_[i]->clear();
Francois Gygi committed
241

242
  int size=amat_[0]->size();
Francois Gygi committed
243

244 245
  // allocate matrix for products of projected parts
  DoubleMatrix products(c.context(),c.n(),c.n(),c.nb(),c.nb());
Francois Gygi committed
246

247 248 249 250 251 252
  for ( int i_proj=0; i_proj<rmat_.size(); i_proj++ )
  {
    // get index of the subbox
    int i_slice[2];
    i_slice[1]= i_proj/ndiv_[0];
    i_slice[0]= i_proj-ndiv_[0]*i_slice[1];
Francois Gygi committed
253

254 255
    // compute product of projections in real space
    products.gemm('t','n',1.0,(*rmat_[i_proj]),(*rmat_[i_proj]),0.0);
Francois Gygi committed
256

257 258 259 260 261 262 263
    // add product to the x/y A matrices
    for ( int ilevel = 0, imat = 0; ilevel < nlevelsmax_; ilevel++ )
    {
      for ( int idir = 0; idir < 3; idir++ )
      {
        // matrix A is an x or y matrix
        if ( ilevel<nlevels_[idir] && idir<2 )
Francois Gygi committed
264
        {
265
          if ( walsh( ilevel , ndiv_[idir] , i_slice[idir] ) )
Francois Gygi committed
266
          {
267 268 269 270 271
            // add product to the matrix
            double *coeff_source=products.valptr(0);
            double *coeff_destination=amat_[imat]->valptr(0);
            for ( int i=0; i<size; i++ )
              coeff_destination[i]+=coeff_source[i];
Francois Gygi committed
272
          }
273 274 275 276 277 278
          imat++;
        }
        // else: matrix A is a z matrix
        else if ( ilevel<nlevels_[idir] )
        {
          imat++;
Francois Gygi committed
279 280
        }
      }
281 282
    }
  }
Francois Gygi committed
283

284 285 286 287 288 289 290
  // normalize xy matrices
  for ( int ilevel = 0, imat = 0; ilevel<nlevelsmax_; ilevel++ )
  {
    for ( int idir = 0; idir < 3; idir++ )
    {
      // matrix A is a projector in the  x or y direction
      if ( ilevel < nlevels_[idir] && idir<2 )
Francois Gygi committed
291
      {
292 293 294 295 296 297 298 299 300 301 302
        // normalize coeffs
        double *coeff=amat_[imat]->valptr(0);
        double fac = 1.0 / (np_[0]*np_[1]*np_[2]);
        for ( int i = 0; i < size; i++ )
          coeff[i] *= fac;
        imat++;
      }
      // else: matrix A is a projector in the z direction
      else if ( ilevel<nlevels_[idir] )
      {
        imat++;
Francois Gygi committed
303 304
      }
    }
305 306
  }

Francois Gygi committed
307
#ifdef TIMING
308
  tmap["xy products"].stop();
Francois Gygi committed
309 310
#endif

311
  // compute the z A matrices
Francois Gygi committed
312
#ifdef TIMING
313
  tmap["z products"].start();
Francois Gygi committed
314
#endif
315 316 317 318 319 320 321 322 323 324 325
  int imat=0;
  // matrix for projected wf
  ComplexMatrix c_pz(c.context(),c.m(),c.n(),c.mb(),c.nb());;
  // proxy for complex->real matrix product
  DoubleMatrix c_proxy(c);
  DoubleMatrix c_pz_proxy(c_pz);

  for ( int ilevel=0; ilevel<nlevelsmax_; ilevel++ )
  {
    // loop on directions
    for ( int idir=0; idir<3; idir++ )
Francois Gygi committed
326
    {
327 328
      // matrix A is for x or y direction
      if ( ilevel<nlevels_[idir] && idir<2 )
Francois Gygi committed
329
      {
330 331 332 333 334 335
        imat++;
      }
      // else: matrix A is for z direction
      else if ( ilevel<nlevels_[idir] )
      {
        for ( int n = 0; n < nstloc_; n++ )
Francois Gygi committed
336
        {
337 338 339
          // p_rmat: pointers to rmat arrays
          vector<double *> p_rmat( ndiv_[0]*ndiv_[1] );
          for ( int iproj=0; iproj<ndiv_[0]*ndiv_[1]; iproj++ )
Francois Gygi committed
340
          {
341 342
            int index = n * rmat_[iproj]->mloc();
            p_rmat[iproj] = rmat_[iproj]->valptr(index);
Francois Gygi committed
343
          }
344 345
          // save values of wf*walsh in rmat arrays
          for ( int iz = 0; iz < np2loc_; iz++ )
Francois Gygi committed
346
          {
347 348 349
            int izglobal = ft_->np2_first() + iz;
            double walsh_z = walsh(ilevel, np_[2], izglobal);
            for ( int ixy = 0; ixy < np01_; ixy++ )
Francois Gygi committed
350
            {
351 352 353 354
              int i = ixy + iz * np01_;
              int xy_proj = xy_proj_[ixy];
              wftmp[i] = (*p_rmat[xy_proj]) * walsh_z;
              p_rmat[xy_proj]++;
Francois Gygi committed
355 356
            }
          }
357
          ft_->forward(&wftmp[0],c_pz.valptr(c_pz.mloc()*n));
Francois Gygi committed
358
        }
359 360 361 362 363 364
        // compute the product
        // factor -2.0 in next line: G and -G
        amat_[imat]->gemm('t','n',2.0,c_pz_proxy,c_proxy,0.0);
        // rank-1 update using first row of cd_proxy() and c_proxy
        amat_[imat]->ger(-1.0,c_pz_proxy,0,c_proxy,0);
        imat++;
Francois Gygi committed
365 366
      }
    }
367
  }
Francois Gygi committed
368
#ifdef TIMING
369
  tmap["z products"].stop();
Francois Gygi committed
370 371 372
#endif

#ifdef DEBUG
373
  // check the values of the amat matrices
Francois Gygi committed
374
#ifdef TIMING
375
  tmap["check_amat"].start();
Francois Gygi committed
376
#endif
377
  check_amat(c);
Francois Gygi committed
378
#ifdef TIMING
379
  tmap["check_amat"].stop();
Francois Gygi committed
380 381 382
#endif
#endif

383 384 385
  // set to zero matrix elements of the matrices amat_[i] if they couple
  // states with differing occupation numbers
  trim_amat(occ);
Francois Gygi committed
386 387

#ifdef DEBUG_PRINT_MAT
388 389 390 391 392
  for ( int k = 0; k < amat_.size(); k++ )
  {
    cout << "A(k=" << k << "):" << endl;
    cout << *amat_[k];
  }
Francois Gygi committed
393 394
#endif

395
  // joint approximate diagonalization of the matrices amat (jade)
Francois Gygi committed
396
#ifdef TIMING
397
  tmap["jade"].start();
Francois Gygi committed
398 399
#endif

400 401
  // diagonal values adiag_[k][i]
  // adiag_ is resized by jade
Francois Gygi committed
402

403
  // diagonalize projectors
Francois Gygi committed
404 405
  // int nsweep = jade(maxsweep,tol,amat_,*u_,adiag_);
  jade(maxsweep,tol,amat_,*u_,adiag_);
406
#ifdef TIMING
407
  if ( ctxt_.onpe0() )
408 409 410
    cout << "Bisection::compute_transform: nsweep=" << nsweep
         << " maxsweep=" << maxsweep << " tol=" << tol << endl;
#endif
Francois Gygi committed
411 412

#ifdef TIMING
413
  tmap["jade"].stop();
Francois Gygi committed
414 415
#endif

416 417 418
#ifdef DEBUG_PRINT_MAT
  cout << "U:" << endl;
  cout << *u_;
Francois Gygi committed
419 420
#endif

421 422 423 424
  for ( int imat=0; imat<nmat_; imat++ )
  {
    // broadcast diagonal of all matrices a to all tasks
    MPI_Bcast( (void *) &adiag_[imat][0], adiag_[imat].size(),
425
               MPI_DOUBLE, 0, ctxt_.comm() );
426
  }
Francois Gygi committed
427 428 429 430 431 432 433
  // print timers
#ifdef TIMING
  for ( map<string,Timer>::iterator i = tmap.begin(); i != tmap.end(); i++ )
  {
    double time = (*i).second.real();
    double tmin = time;
    double tmax = time;
434 435 436
    ctxt_.dmin(1,1,&tmin,1);
    ctxt_.dmax(1,1,&tmax,1);
    if ( ctxt_.myproc()==0 )
Francois Gygi committed
437
    {
438 439
      string s = "name=\"" + (*i).first + "\"";
      cout << "<timing " << left << setw(22) << s
440 441
           << " min=\"" << setprecision(3) << tmin << "\""
           << " max=\"" << setprecision(3) << tmax << "\"/>"
Francois Gygi committed
442 443 444 445 446 447 448
           << endl;
    }
  }
#endif
}

////////////////////////////////////////////////////////////////////////////////
449
void Bisection::compute_localization(double epsilon)
Francois Gygi committed
450
{
451 452
  // compute localization vector for a threshold epsilon
  for ( int n = 0; n < nst_; n++ )
Francois Gygi committed
453
  {
454 455
    localization_[n] = 0;
    for ( int imat = 0; imat < nmat_; imat++ )
Francois Gygi committed
456
    {
457 458 459 460 461 462
      if ( adiag_[imat][n] < epsilon )
        localization_[n] += 1<<(2*imat);
      else if ( adiag_[imat][n] > 1.0-epsilon)
        localization_[n] += 1<<(2*imat+1);
      else
        localization_[n] += (1<<(2*imat)) + (1<<(2*imat+1));
Francois Gygi committed
463 464 465
    }
  }

466 467 468
#ifdef DEBUG
  // print localization vector and number of overlaps (including self)
  // for each state
469
  if ( ctxt_.onpe0() )
Francois Gygi committed
470
  {
471 472
    int sum = 0;
    for ( int i = 0; i < nst_; i++ )
Francois Gygi committed
473
    {
474 475 476 477 478 479 480 481 482 483 484
      int count = 0;
      for ( int j = 0; j < nst_; j++ )
      {
        if ( overlap(i,j) )
          count++;
      }
      cout << "localization[" << i << "]: "
           << localization_[i] << " "
           << bitset<30>(localization_[i]) << "  overlaps: "
           << count << endl;
      sum += count;
Francois Gygi committed
485
    }
486 487
    cout << "total overlaps: " << sum << " / " << nst_*nst_
         << " = " << ((double) sum)/(nst_*nst_) << endl;
Francois Gygi committed
488
  }
489
#endif
Francois Gygi committed
490

491 492
  // broadcast localization to all tasks to ensure consistency
  MPI_Bcast( (void *) &localization_[0], localization_.size(),
493
             MPI_LONG, 0, ctxt_.comm() );
Francois Gygi committed
494 495 496 497

}

////////////////////////////////////////////////////////////////////////////////
498
void Bisection::forward(SlaterDet& sd)
Francois Gygi committed
499
{
500 501
  forward(*u_,sd);
}
Francois Gygi committed
502

503 504 505 506 507 508 509 510 511 512 513
////////////////////////////////////////////////////////////////////////////////
void Bisection::forward(DoubleMatrix& u, SlaterDet& sd)
{
  // apply the bisection transformation to the SlaterDet sd
  // apply the rotation u to sd.c()
  ComplexMatrix& c = sd.c();
  ComplexMatrix cp(c);
  DoubleMatrix cp_proxy(cp);
  DoubleMatrix c_proxy(c);
  c_proxy.gemm('n','n',1.0,cp_proxy,u,0.0);
}
Francois Gygi committed
514

515 516 517 518 519
////////////////////////////////////////////////////////////////////////////////
void Bisection::backward(SlaterDet& sd)
{
  backward(*u_,sd);
}
Francois Gygi committed
520

521 522 523 524 525 526 527 528 529 530
////////////////////////////////////////////////////////////////////////////////
void Bisection::backward(DoubleMatrix& u, SlaterDet& sd)
{
  // apply the inverse bisection transformation to SlaterDet sd
  // apply rotation u^T to sd
  ComplexMatrix& c = sd.c();
  ComplexMatrix cp(c);
  DoubleMatrix cp_proxy(cp);
  DoubleMatrix c_proxy(c);
  c_proxy.gemm('n','t',1.0,cp_proxy,u,0.0);
Francois Gygi committed
531 532 533
}

////////////////////////////////////////////////////////////////////////////////
534
bool Bisection::check_amat(const ComplexMatrix &c)
Francois Gygi committed
535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677
{
  // allocate memory for copy of the wave functions
  ComplexMatrix cd(c.context(),c.m(),c.n(),c.mb(),c.nb());

  // precompute values of Walsh functions in three directions
  // wx[l][i] = value of level l Walsh function at position i in direction x
  vector<vector<vector<int> > > w(3);
  for ( int idir=0; idir<3; idir++ )
  {
    w[idir].resize(nlevels_[idir]);
    //
    for ( int l = 0; l < nlevels_[idir]; l++ )
    {
      w[idir][l].resize(np_[idir]);
      for ( int i = 0; i < np_[idir]; i++ )
      {
        w[idir][l][i] = walsh(l,np_[idir],i);
      }
    }
  }

  // compute matrices B_k = <wf|dwf>
  vector<DoubleMatrix*> bmat(nmat_);
  for ( int k = 0; k < bmat.size(); k++ )
  {
    bmat[k] = new DoubleMatrix(c.context(),c.n(),c.n(),c.nb(),c.nb());
  }
  DoubleMatrix c_proxy(c);
  DoubleMatrix cd_proxy(cd);

  vector<complex<double> > wftmp(ft_->np012loc());
  complex<double> *f = &wftmp[0];

  // compute matrices A at all levels
  for ( int l=0 , imat=0; l < nlevelsmax_; l++ )
  {
    // x direction
    if ( l<nlevels_[0] )
    {
      cd_proxy = c_proxy;
      for ( int n = 0; n < c.nloc(); n++ )
      {
        for ( int i = 0; i < np012loc_; i++ )
          f[i] = 0.0;
        ft_->backward(cd.cvalptr(cd.mloc()*n),&wftmp[0]);
        for ( int i = 0; i < np_[0]; i++ )
        {
          if ( w[0][l][i] == 0 )
          {
            for ( int j = 0; j < np_[1]; j++ )
              for ( int k = 0; k < ft_->np2_loc(); k++ )
                f[i +  np_[0] * ( j + np_[1] * k )] = 0.0;
          }
        }
        ft_->forward(&wftmp[0],cd.valptr(cd.mloc()*n));
      }
      // factor -2.0 in next line: G and -G
      bmat[imat]->gemm('t','n',2.0,cd_proxy,c_proxy,0.0);
      // rank-1 update using first row of cd_proxy() and c_proxy
      bmat[imat]->ger(-1.0,cd_proxy,0,c_proxy,0);
      imat++;
    }

    // y direction
    if ( l<nlevels_[1] )
    {
      cd_proxy = c_proxy;
      for ( int n = 0; n < c.nloc(); n++ )
      {
        for ( int i = 0; i < np012loc_; i++ )
          f[i] = 0.0;
        ft_->backward(cd.cvalptr(cd.mloc()*n),&wftmp[0]);
        for ( int j = 0; j < np_[1]; j++ )
        {
          if ( w[1][l][j] == 0 )
          {
            for ( int i = 0; i < np_[0]; i++ )
              for ( int k = 0; k < ft_->np2_loc(); k++ )
                f[i +  np_[0] * ( j + np_[1] * k )] = 0.0;
          }
        }
        ft_->forward(&wftmp[0],cd.valptr(cd.mloc()*n));
      }
      // factor -2.0 in next line: G and -G
      bmat[imat]->gemm('t','n',2.0,cd_proxy,c_proxy,0.0);
      // rank-1 update using first row of cd_proxy() and c_proxy
      bmat[imat]->ger(-1.0,cd_proxy,0,c_proxy,0);
      imat++;
    }

    // z direction
    if ( l<nlevels_[2] )
    {
      cd_proxy = c_proxy;
      for ( int n = 0; n < c.nloc(); n++ )
      {
        for ( int i = 0; i < np012loc_; i++ )
          f[i] = 0.0;
        ft_->backward(cd.cvalptr(cd.mloc()*n),&wftmp[0]);
        for ( int k = 0; k < ft_->np2_loc(); k++ )
        {
          int kglobal = ft_->np2_first() + k;
          const int istart = k * np_[0] * np_[1];
          if ( w[2][l][kglobal] == 0 )
          {
            for ( int ij = 0; ij < np_[0]*np_[1]; ij++ )
              f[istart+ij] = 0.0;
          }
        }
        ft_->forward(&wftmp[0],cd.valptr(cd.mloc()*n));
      }
      // factor -2.0 in next line: G and -G
      bmat[imat]->gemm('t','n',2.0,cd_proxy,c_proxy,0.0);
      // rank-1 update using first row of cd_proxy() and c_proxy
      bmat[imat]->ger(-1.0,cd_proxy,0,c_proxy,0);
      imat++;
    }
  } // for l

  // testing the matrices
  for ( int imat=0; imat<nmat_; imat++ )
  {
    double *a=amat_[imat]->valptr(0);
    double *b=bmat[imat]->valptr(0);
    int ncoeff=amat_[imat]->mloc()*amat_[imat]->nloc();

    for ( int i=0; i<ncoeff; i++ )
    {
      if (fabs(a[i]-b[i])>1e-10)
      {
        cout << "error > 1.e-10 for matrix " << imat << endl;
      }
    }
  }
  return true;
}

////////////////////////////////////////////////////////////////////////////////
void Bisection::trim_amat(const vector<double>& occ)
{
  // set to zero the matrix elements of the matrices amat_[k] if they couple
  // states with differing occupation numbers

Francois Gygi committed
678
  const double trim_tol = 1.e-6;
Francois Gygi committed
679 680 681 682 683 684 685 686
  // check if all occupation numbers are the same
  double occ_max = 0.0, occ_min = 2.0;
  for ( int i = 0; i < occ.size(); i++ )
  {
    occ_max = max(occ_max, occ[i]);
    occ_min = min(occ_min, occ[i]);
  }
  // return if all occupation numbers are equal
Francois Gygi committed
687
  if ( fabs(occ_max-occ_min) < trim_tol )
Francois Gygi committed
688 689 690 691 692 693 694 695 696 697 698 699 700
    return;

  const int mloc = amat_[0]->mloc();
  const int nloc = amat_[0]->nloc();
  // loop over elements local to this task
  for ( int i = 0; i < mloc; i++ )
  {
    const int iglobal = amat_[0]->iglobal(i);
    for ( int j = 0; j < nloc; j++ )
    {
      const int jglobal = amat_[0]->jglobal(j);

      const int ival = i + mloc * j;
Francois Gygi committed
701
      if ( fabs(occ[iglobal] - occ[jglobal]) > trim_tol )
Francois Gygi committed
702 703 704 705 706 707 708
        for ( int k = 0; k < amat_.size(); k++ )
          (*amat_[k])[ival] = 0.0;
    }
  }
}

////////////////////////////////////////////////////////////////////////////////
709
bool Bisection::overlap(int i, int j) const
710 711 712 713 714
{
  return overlap(localization_,i,j);
}

////////////////////////////////////////////////////////////////////////////////
715
bool Bisection::overlap(const vector<long int>& loc_, int i, int j) const
Francois Gygi committed
716
{
717 718 719 720 721
  // overlap: return true if the functions i and j overlap according
  // to the localization vector loc_
  long int loc_i = loc_[i];
  long int loc_j = loc_[j];
  while ( loc_i!=0 && loc_j!=0 )
Francois Gygi committed
722 723
  {
    // get the weight of projections for each state
724 725 726 727
    bool p_right_i = loc_i & 1;
    bool p_left_i  = loc_i & 2;
    bool p_right_j = loc_j & 1;
    bool p_left_j  = loc_j & 2;
Francois Gygi committed
728 729 730 731 732

    // return false as soon as the states are found to be separated
    if ( !( ( p_right_i && p_right_j ) || ( p_left_i && p_left_j ) ) )
      return false;

733 734
    loc_i >>= 2;
    loc_j >>= 2;
Francois Gygi committed
735 736 737 738 739 740 741
  }

  // return true if the states overlap
  return true;
}

////////////////////////////////////////////////////////////////////////////////
742
double Bisection::pair_fraction(void) const
Francois Gygi committed
743 744
{
  // pair_fraction: return fraction of pairs having non-zero overlap
745
  // count pairs (i,j) having non-zero overlap for i != j only
Francois Gygi committed
746
  int sum = 0;
747
  for ( int i = 1; i < nst_; i++ )
Francois Gygi committed
748 749
  {
    int count = 0;
750
    for ( int j = i+1; j < nst_; j++ )
Francois Gygi committed
751
    {
752
      if ( overlap(i,j) )
Francois Gygi committed
753 754 755 756
        count++;
    }
    sum += count;
  }
757 758 759
  // add overlap with self: (i,i)
  sum += nst_;
  return ((double) sum)/((nst_*(nst_+1))/2);
Francois Gygi committed
760 761 762
}

////////////////////////////////////////////////////////////////////////////////
763
double Bisection::size(int i) const
Francois Gygi committed
764 765
{
  // size: return fraction of the domain on which state i is non-zero
766
  long int loc = localization_[i];
Francois Gygi committed
767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785
  double size = 1.0;
  // process all projectors
  while ( loc != 0 )
  {
    // weight of projections
    bool p_right = loc & 1;
    bool p_left  = loc & 2;

    if ( (p_right && !p_left) || (!p_right && p_left) )
      size *= 0.5;

    loc >>= 2;
  }

  // return true if the states overlap
  return size;
}

////////////////////////////////////////////////////////////////////////////////
786
double Bisection::total_size(void) const
Francois Gygi committed
787 788 789
{
  // total_size: return normalized sum of sizes
  double sum = 0.0;
790 791 792
  for ( int i = 0; i < nst_; i++ )
    sum += size(i);
  return sum / nst_;
Francois Gygi committed
793
}