SlaterDet.C 54.7 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 18 19 20 21 22
// SlaterDet.C
//
////////////////////////////////////////////////////////////////////////////////

#include "SlaterDet.h"
#include "FourierTransform.h"
#include "Context.h"
#include "blas.h" // daxpy
Francois Gygi committed
23
#include "Base64Transcoder.h"
24
#include "SharedFilePtr.h"
Francois Gygi committed
25
#include "Timer.h"
Francois Gygi committed
26 27

#include <cstdlib>
28
#include <cstring> // memcpy
Francois Gygi committed
29 30 31
#include <iostream>
#include <iomanip>
#include <sstream>
Francois Gygi committed
32
#include <limits>
Francois Gygi committed
33 34 35
using namespace std;

////////////////////////////////////////////////////////////////////////////////
36
SlaterDet::SlaterDet(const Context& ctxt, D3vector kpoint) : ctxt_(ctxt),
37
 c_(ctxt)
Francois Gygi committed
38
{
39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56
  // create cartesian communicator mapped on ctxt
  int ndims=2;
  // Note: MPI_Cart_comm uses row-major ordering. Need to give
  // transposed dimensions as input arguments
  int dims[2] = {ctxt.npcol(), ctxt.nprow()};
  int periods[2] = { 0, 0};
  int reorder = 0;
  MPI_Comm comm;
  MPI_Cart_create(ctxt.comm(),ndims,dims,periods,reorder,&comm);

  int size, myrank;
  MPI_Comm_size(comm,&size);
  MPI_Comm_rank(comm,&myrank);
  int coords[2];
  MPI_Cart_coords(comm,myrank,2,coords);

#ifdef DEBUG
  for ( int i = 0; i < size; i++ )
57
  {
58 59 60 61
    MPI_Barrier(comm);
    if ( myrank == i )
      cout << myrank << ": myrow=" << ctxt.myrow() << " mycol=" << ctxt.mycol()
           << " coords= " << coords[0] << ", " << coords[1] << endl;
62
  }
63 64 65 66 67 68 69 70 71 72 73 74 75
#endif

  // Split the cartesian communicator comm to define my_col_comm_
  // MPI_Cart_create uses row-major ordering. Need to keep the second
  // dimension to get a communicator corresponding to a column of ctxt
  int remain_dims[2] = { 0, 1 };
  MPI_Cart_sub(comm, remain_dims, &my_col_comm_);

  // Free the cartesian communicator
  MPI_Comm_free(&comm);

  // define basis on the column subcommunicator
  basis_ = new Basis(my_col_comm_,kpoint);
Francois Gygi committed
76
}
Francois Gygi committed
77 78 79

////////////////////////////////////////////////////////////////////////////////
SlaterDet::SlaterDet(const SlaterDet& rhs) : ctxt_(rhs.context()),
80 81 82 83
  basis_(new Basis(*(rhs.basis_))), c_(rhs.c_)
  {
    MPI_Comm_dup(rhs.my_col_comm_,&my_col_comm_);
  }
84

Francois Gygi committed
85
////////////////////////////////////////////////////////////////////////////////
86
SlaterDet::~SlaterDet()
87 88
{
  delete basis_;
Francois Gygi committed
89
  // cout << ctxt_.mype() << ": SlaterDet::dtor: ctxt=" << ctxt_;
90 91 92 93 94 95 96 97 98 99
#ifdef TIMING
  for ( TimerMap::iterator i = tmap.begin(); i != tmap.end(); i++ )
  {
    double time = (*i).second.real();
    double tmin = time;
    double tmax = time;
    ctxt_.dmin(1,1,&tmin,1);
    ctxt_.dmax(1,1,&tmax,1);
    if ( ctxt_.myproc()==0 )
    {
100 101
      string s = "name=\"" + (*i).first + "\"";
      cout << "<timing " << left << setw(22) << s
102 103
           << " min=\"" << setprecision(3) << tmin << "\""
           << " max=\"" << setprecision(3) << tmax << "\"/>"
104 105 106 107
           << endl;
    }
  }
#endif
108
  MPI_Comm_free(&my_col_comm_);
Francois Gygi committed
109 110 111
}

////////////////////////////////////////////////////////////////////////////////
112
void SlaterDet::resize(const UnitCell& cell, const UnitCell& refcell,
Francois Gygi committed
113 114
  double ecut, int nst)
{
115 116
  // Test in next line should be replaced by test on basis min/max indices
  // to signal change in basis vectors
117
  //if ( basis_->refcell().volume() != 0.0 && !refcell.encloses(cell) )
118 119
  //{
    //cout << " SlaterDet::resize: cell=" << cell;
120
    //cout << " SlaterDet::resize: refcell=" << basis_->refcell();
121
    //throw SlaterDetException("could not resize: cell not in refcell");
122
  //}
123

Francois Gygi committed
124 125
  try
  {
126
    // create a temporary copy of the basis
127 128 129
    Basis btmp(*basis_);

    // perform normal resize operations, possibly resetting contents of c_
130
    basis_->resize(cell,refcell,ecut);
Francois Gygi committed
131 132
    occ_.resize(nst);
    eig_.resize(nst);
133

134
    const int mb = basis_->maxlocalsize();
135 136
    const int m = ctxt_.nprow() * mb;
    const int nb = nst/ctxt_.npcol() + (nst%ctxt_.npcol() > 0 ? 1 : 0);
137

138 139 140 141 142 143
    // check if the dimensions of c_ must change
    if ( m!=c_.m() || nst!=c_.n() || mb!=c_.mb() || nb!=c_.nb() )
    {
      // make a copy of c_ before resize
      ComplexMatrix ctmp(c_);
      c_.resize(m,nst,mb,nb);
Francois Gygi committed
144
      init();
145

146 147 148
      // check if data can be copied from temporary copy
      // It is assumed that nst and ecut are not changing at the same time
      // Only the cases where one change at a time occurs is covered
149

150 151
      // consider only cases where the dimensions are finite
      if ( c_.m() > 0 && c_.n() > 0 )
152
      {
153 154
        // first case: only nst has changed
        if ( c_.m() == ctmp.m() && c_.n() != ctmp.n() )
155
        {
156 157 158 159 160 161 162 163 164 165 166 167 168 169 170
          //cout << "SlaterDet::resize: c_m/n=   "
          //     << c_.m() << "/" << c_.n() << endl;
          //cout << "SlaterDet::resize: ctmp_m/n=" << ctmp.m()
          //     << "/" << ctmp.n() << endl;
          // nst has changed, basis is unchanged
          // copy coefficients up to min(n_old, n_new)
          if ( c_.n() < ctmp.n() )
          {
            c_.getsub(ctmp,ctmp.m(),c_.n(),0,0);
          }
          else
          {
            c_.getsub(ctmp,ctmp.m(),ctmp.n(),0,0);
          }
          gram();
171
        }
172 173 174
        // second case: basis was resized, nst unchanged
        if ( btmp.ecut() > 0.0 && basis_->ecut() > 0.0 &&
             c_.m() != ctmp.m() && c_.n() == ctmp.n() )
175
        {
176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195
          // transform all states to real space and interpolate
          int np0 = max(basis_->np(0),btmp.np(0));
          int np1 = max(basis_->np(1),btmp.np(1));
          int np2 = max(basis_->np(2),btmp.np(2));
          //cout << " SlaterDet::resize: grid: np0/1/2: "
          //     << np0 << " " << np1 << " " << np2 << endl;
          // FourierTransform tf1(oldbasis, new basis grid)
          // FourierTransform tf2(newbasis, new basis grid)
          FourierTransform ft1(btmp,np0,np1,np2);
          FourierTransform ft2(*basis_,np0,np1,np2);
          // allocate real-space grid
          valarray<complex<double> > tmpr(ft1.np012loc());
          // transform each state from old basis to grid to new basis
          for ( int n = 0; n < nstloc(); n++ )
          {
            for ( int i = 0; i < tmpr.size(); i++ )
              tmpr[i] = 0.0;
            ft1.backward(ctmp.cvalptr(n*ctmp.mloc()),&tmpr[0]);
            ft2.forward(&tmpr[0], c_.valptr(n*c_.mloc()));
          }
196 197 198
        }
      }
    }
Francois Gygi committed
199 200 201 202 203 204
  }
  catch ( bad_alloc )
  {
    cout << " bad_alloc exception caught in SlaterDet::resize" << endl;
    throw;
  }
Francois Gygi committed
205
#if 0
206 207 208 209 210 211 212 213 214
  // print error in imaginary part of c(G=0)
  double imag_err = g0_imag_error();
  if ( ctxt_.onpe0() )
  {
    cout.setf(ios::scientific,ios::floatfield);
    cout << " SlaterDet::resize: imag error on exit: " << imag_err << endl;
  }
#endif
  cleanup();
Francois Gygi committed
215
}
216
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
217
void SlaterDet::init(void)
218 219
{
  // initialize coefficients with lowest plane waves
220
  if ( c_.n() <= basis_->size() )
221 222 223 224
  {
    // initialize c_
    c_.clear();
    const double s2i = 1.0 / sqrt(2.0);
225

226 227 228 229 230 231 232 233 234
    // for each n, find the smallest g vector and initialize
    int ismallest = 0;
    // on each process, basis.isort(ismallest) is the index of the smallest
    // local g vector
    for ( int n = 0; n < c_.n(); n++ )
    {
      double value = 1.0;
      if ( basis().real() && n != 0 )
        value = s2i;
235

236
      // find process row holding the smallest g vector
237 238 239
      // kpg2: size^2 of smallest vector on this task
      // set kpg2 to largest double value if localsize == 0
      double kpg2 = numeric_limits<double>::max();
240
      if ( ismallest < basis_->localsize() )
241 242 243
      {
        kpg2 = basis_->kpg2(basis_->isort(ismallest));
      }
Francois Gygi committed
244
      // cout << "smallest vector on proc " << ctxt_.mype()
245
      //      << " has norm " << kpg2 << endl;
246
      int minrow, mincol;
247 248
      ctxt_.dmin('c',' ',1,1,&kpg2,1,&minrow,&mincol,1,-1,-1);

249 250 251 252 253
      // find column hosting state n
      int pc = c_.pc(n);
      int pr = minrow;
      if ( pr == ctxt_.myrow() )
      {
254
        int iii = basis_->isort(ismallest);
255 256 257 258 259
        ismallest++; // increment on entire process row
        if ( pc == ctxt_.mycol() )
        {
          // cout << " n=" << n << " on process "
          //      << pr << "," << pc
260 261 262
          //      << " vector " << basis_->idx(3*iii) << " "
          //      << basis_->idx(3*iii+1) << " "
          //      << basis_->idx(3*iii+2) << " norm="
Francois Gygi committed
263
          //      << basis_->g2(iii) << " "
264 265 266 267 268 269 270 271 272
          //      << value << endl;
          int jjj = c_.m(n) * c_.nb() + c_.y(n);
          int index = iii+c_.mloc()*jjj;
          c_[index] = complex<double> (value,0.0);
        }
      }
    }
  }
}
Francois Gygi committed
273 274

////////////////////////////////////////////////////////////////////////////////
275
void SlaterDet::compute_density(FourierTransform& ft,
Francois Gygi committed
276 277
  double weight, double* rho) const
{
Francois Gygi committed
278
  //Timer tm_ft, tm_rhosum;
Francois Gygi committed
279 280 281
  // compute density of the states residing on my column of ctxt_
  assert(occ_.size() == c_.n());
  vector<complex<double> > tmp(ft.np012loc());
282

283 284
  assert(basis_->cell().volume() > 0.0);
  const double omega_inv = 1.0 / basis_->cell().volume();
Francois Gygi committed
285
  const int np012loc = ft.np012loc();
286

287
  if ( basis_->real() )
Francois Gygi committed
288 289 290 291 292 293
  {
    // transform two states at a time
    for ( int n = 0; n < nstloc()-1; n++, n++ )
    {
      // global n index
      const int nn = ctxt_.mycol() * c_.nb() + n;
Francois Gygi committed
294 295
      const double fac1 = weight * omega_inv * occ_[nn];
      const double fac2 = weight * omega_inv * occ_[nn+1];
296

Francois Gygi committed
297 298
      if ( fac1 + fac2 > 0.0 )
      {
Francois Gygi committed
299
        //tm_ft.start();
Francois Gygi committed
300
        ft.backward(c_.cvalptr(n*c_.mloc()),
Francois Gygi committed
301 302
                    c_.cvalptr((n+1)*c_.mloc()),&tmp[0]);
        //tm_ft.stop();
Francois Gygi committed
303 304
        const double* psi = (double*) &tmp[0];
        int ii = 0;
Francois Gygi committed
305
        //tm_rhosum.start();
Francois Gygi committed
306 307 308 309 310 311 312
        for ( int i = 0; i < np012loc; i++ )
        {
          const double psi1 = psi[ii];
          const double psi2 = psi[ii+1];
          rho[i] += fac1 * psi1 * psi1 + fac2 * psi2 * psi2;
          ii++; ii++;
        }
Francois Gygi committed
313
        //tm_rhosum.start();
Francois Gygi committed
314 315 316 317 318 319 320
      }
    }
    if ( nstloc() % 2 != 0 )
    {
      const int n = nstloc()-1;
      // global n index
      const int nn = ctxt_.mycol() * c_.nb() + n;
Francois Gygi committed
321
      const double fac1 = weight * omega_inv * occ_[nn];
322

Francois Gygi committed
323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343
      if ( fac1 > 0.0 )
      {
        ft.backward(c_.cvalptr(n*c_.mloc()),&tmp[0]);
        const double* psi = (double*) &tmp[0];
        int ii = 0;
        for ( int i = 0; i < np012loc; i++ )
        {
          const double psi1 = psi[ii];
          rho[i] += fac1 * psi1 * psi1;
          ii++; ii++;
        }
      }
    }
  }
  else
  {
    // only one transform at a time
    for ( int n = 0; n < nstloc(); n++ )
    {
      // global n index
      const int nn = ctxt_.mycol() * c_.nb() + n;
Francois Gygi committed
344
      const double fac = weight * omega_inv * occ_[nn];
345

Francois Gygi committed
346 347 348 349 350 351 352 353
      if ( fac > 0.0 )
      {
        ft.backward(c_.cvalptr(n*c_.mloc()),&tmp[0]);
        for ( int i = 0; i < np012loc; i++ )
          rho[i] += fac * norm(tmp[i]);
      }
    }
  }
354
  // cout << "SlaterDet: compute_density: ft_bwd time: "
Francois Gygi committed
355
  //      << tm_ft.real() << endl;
356
  // cout << "SlaterDet: compute_density: rhosum time: "
Francois Gygi committed
357
  //      << tm_rhosum.real() << endl;
Francois Gygi committed
358 359 360
}

////////////////////////////////////////////////////////////////////////////////
361
void SlaterDet::rs_mul_add(FourierTransform& ft,
Francois Gygi committed
362
  const double* v, SlaterDet& sdp) const
Francois Gygi committed
363 364 365 366
{
  // transform states to real space, multiply states by v[r] in real space
  // transform back to reciprocal space and add to sdp
  // sdp[n] += v * sd[n]
367

Francois Gygi committed
368 369
  vector<complex<double> > tmp(ft.np012loc());
  vector<complex<double> > ctmp(2*c_.mloc());
370

Francois Gygi committed
371 372 373 374
  const int np012loc = ft.np012loc();
  const int mloc = c_.mloc();
  double* dcp = (double*) sdp.c().valptr();

375
  if ( basis_->real() )
Francois Gygi committed
376 377 378 379 380 381
  {
    // transform two states at a time
    for ( int n = 0; n < nstloc()-1; n++, n++ )
    {
      ft.backward(c_.cvalptr(n*mloc),
                 c_.cvalptr((n+1)*mloc),&tmp[0]);
Francois Gygi committed
382 383

      #pragma omp parallel for
Francois Gygi committed
384
      for ( int i = 0; i < np012loc; i++ )
Francois Gygi committed
385 386
        tmp[i] *= v[i];

Francois Gygi committed
387 388 389 390
      ft.forward(&tmp[0], &ctmp[0], &ctmp[mloc]);
      int len = 4 * mloc;
      int inc1 = 1;
      double alpha = 1.0;
Francois Gygi committed
391
      daxpy(&len,&alpha,(double*)&ctmp[0],&inc1,&dcp[2*n*mloc],&inc1);
Francois Gygi committed
392 393 394 395 396
    }
    if ( nstloc() % 2 != 0 )
    {
      const int n = nstloc()-1;
      ft.backward(c_.cvalptr(n*mloc),&tmp[0]);
Francois Gygi committed
397 398

      #pragma omp parallel for
Francois Gygi committed
399
      for ( int i = 0; i < np012loc; i++ )
Francois Gygi committed
400 401
        tmp[i] *= v[i];

Francois Gygi committed
402 403 404 405
      ft.forward(&tmp[0], &ctmp[0]);
      int len = 2 * mloc;
      int inc1 = 1;
      double alpha = 1.0;
Francois Gygi committed
406
      daxpy(&len,&alpha,(double*)&ctmp[0],&inc1,&dcp[2*n*mloc],&inc1);
Francois Gygi committed
407 408 409 410 411 412 413 414
    }
  }
  else
  {
    // only one transform at a time
    for ( int n = 0; n < nstloc(); n++ )
    {
      ft.backward(c_.cvalptr(n*mloc),&tmp[0]);
Francois Gygi committed
415 416

      #pragma omp parallel for
Francois Gygi committed
417 418
      for ( int i = 0; i < np012loc; i++ )
        tmp[i] *= v[i];
Francois Gygi committed
419

Francois Gygi committed
420 421 422 423
      ft.forward(&tmp[0], &ctmp[0]);
      int len = 2 * mloc;
      int inc1 = 1;
      double alpha = 1.0;
Francois Gygi committed
424
      daxpy(&len,&alpha,(double*)&ctmp[0],&inc1,&dcp[2*n*mloc],&inc1);
Francois Gygi committed
425 426
    }
  }
427

Francois Gygi committed
428 429 430 431 432
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::gram(void)
{
433
  cleanup();
434
  if ( basis_->real() )
Francois Gygi committed
435 436 437 438 439
  {
    // k = 0 case
    // create a DoubleMatrix proxy for c_
    DoubleMatrix c_proxy(c_);
    DoubleMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
440 441 442
#if TIMING
    tmap["syrk"].start();
#endif
443
    s.syrk('l','t',2.0,c_proxy,0.0);
444 445 446 447
#if TIMING
    tmap["syrk"].stop();
    tmap["syr"].start();
#endif
Francois Gygi committed
448
    s.syr('l',-1.0,c_proxy,0,'r');
449 450 451 452 453
#if TIMING
    tmap["syr"].stop();
    tmap["potrf"].start();
#endif

454 455
#ifdef CHOLESKY_REMAP
    // create a square context for the Cholesky decomposition
456 457
    // int nsq = (int) sqrt((double) ctxt_.size());
    int nsq = CHOLESKY_REMAP;
458 459 460 461 462
    Context csq(nsq,nsq);
    DoubleMatrix ssq(csq,c_.n(),c_.n(),c_.nb(),c_.nb());
    ssq.getsub(s,s.m(),s.n(),0,0);
    ssq.potrf('l'); // Cholesky decomposition: S = L * L^T
    s.getsub(ssq,s.m(),s.n(),0,0);
463 464
#else
    s.potrf('l'); // Cholesky decomposition: S = L * L^T
465
#endif
Francois Gygi committed
466
    // solve triangular system X * L^T = C
467 468 469 470
#if TIMING
    tmap["potrf"].stop();
    tmap["trsm"].start();
#endif
Francois Gygi committed
471
    c_proxy.trsm('r','l','t','n',1.0,s);
472 473 474
#if TIMING
    tmap["trsm"].stop();
#endif
Francois Gygi committed
475 476 477 478 479
  }
  else
  {
    // k != 0 case
    ComplexMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
480
    s.herk('l','c',1.0,c_,0.0);
481 482
    s.potrf('l'); // Cholesky decomposition: S = L * L^H
    // solve triangular system X * L^H = C
Francois Gygi committed
483 484 485 486 487
    c_.trsm('r','l','c','n',1.0,s);
  }
}

////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
488
void SlaterDet::riccati(const SlaterDet& sd)
Francois Gygi committed
489
{
490
  cleanup();
491
  if ( basis_->real() )
Francois Gygi committed
492 493 494 495 496 497
  {
    // k = 0 case
    DoubleMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix r(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    s.identity();
    r.identity();
498

Francois Gygi committed
499 500 501
    DoubleMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix xm(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
502

Francois Gygi committed
503 504 505
    // DoubleMatrix proxy for c_ and sd.c()
    DoubleMatrix c_proxy(c_);
    DoubleMatrix sdc_proxy(sd.c());
506

Francois Gygi committed
507 508 509
#if TIMING
    tmap["riccati_syrk"].start();
#endif
Francois Gygi committed
510 511 512 513
    // Factor -1.0 in next line: -0.5 from definition of s, 2.0 for G and -G
    s.syrk('l','t',-1.0,c_proxy,0.5); // s = 0.5 * ( I - A )
    // symmetric rank-1 update using first row of c_proxy
    s.syr('l',0.5,c_proxy,0,'r');
Francois Gygi committed
514 515 516 517
#if TIMING
    tmap["riccati_syrk"].stop();
    tmap["riccati_gemm"].start();
#endif
Francois Gygi committed
518 519 520 521
    // factor -2.0 in next line: G and -G
    r.gemm('t','n',-2.0,sdc_proxy,c_proxy,1.0); // r = ( I - B )
    // rank-1 update using first row of sdc_proxy() and c_proxy
    r.ger(1.0,sdc_proxy,0,c_proxy,0);
Francois Gygi committed
522 523 524
#if TIMING
    tmap["riccati_gemm"].stop();
#endif
525

Francois Gygi committed
526 527
    xm = s;
    xm.symmetrize('l');
528

Francois Gygi committed
529 530 531
#if TIMING
    tmap["riccati_while"].start();
#endif
Francois Gygi committed
532 533
    s.syrk('l','t',0.5,r,1.0); // s = s + 0.5 * r^T * r
    s.symmetrize('l');
534

Francois Gygi committed
535 536
    double diff = 1.0;
    const double epsilon = 1.e-10;
Francois Gygi committed
537
    const int maxiter = 5;
Francois Gygi committed
538
    int iter = 0;
539

Francois Gygi committed
540 541 542 543 544 545 546
    while ( iter < maxiter && diff > epsilon )
    {
      // x = s - 0.5 * ( r - xm )^T * ( r - xm )
      // Note: t and r are not symmetric, x, xm, and s are symmetric

      for ( int i = 0; i < t.size(); i++ )
        t[i] = r[i] - xm[i];
547

Francois Gygi committed
548 549
      x = s;
      x.syrk('l','t',-0.5,t,1.0);
550

Francois Gygi committed
551 552
      // get full matrix x
      x.symmetrize('l');
553

Francois Gygi committed
554 555
      for ( int i = 0; i < t.size(); i++ )
        t[i] = x[i] - xm[i];
556

Francois Gygi committed
557
      diff = t.nrm2();
558

Francois Gygi committed
559 560 561
      xm = x;
      iter++;
    }
Francois Gygi committed
562 563 564 565
#if TIMING
    tmap["riccati_while"].stop();
    tmap["riccati_symm"].start();
#endif
Francois Gygi committed
566
    c_proxy.symm('r','l',1.0,x,sdc_proxy,1.0);
Francois Gygi committed
567 568 569
#if TIMING
    tmap["riccati_symm"].stop();
#endif
Francois Gygi committed
570 571 572 573 574 575 576 577
  }
  else
  {
    // k != 0 case
    ComplexMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix r(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    s.identity();
    r.identity();
578

Francois Gygi committed
579 580 581
    ComplexMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix xm(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
582

Francois Gygi committed
583 584 585 586
    // s = 0.5 * ( I - A )
    s.herk('l','c',-0.5,c_,0.5);
    // r = ( I - B )
    r.gemm('c','n',-1.0,sd.c(),c_,1.0);
587

Francois Gygi committed
588 589
    xm = s;
    xm.symmetrize('l');
590

Francois Gygi committed
591 592 593
    // s = s + 0.5 * r^H * r
    s.herk('l','c',0.5,r,1.0);
    s.symmetrize('l');
594

Francois Gygi committed
595 596
    double diff = 1.0;
    const double epsilon = 1.e-10;
Francois Gygi committed
597
    const int maxiter = 5;
Francois Gygi committed
598
    int iter = 0;
599

Francois Gygi committed
600 601 602 603 604 605 606
    while ( iter < maxiter && diff > epsilon )
    {
      // x = s - 0.5 * ( r - xm )^H * ( r - xm )
      // Note: t and r are not hermitian, x, xm, and s are hermitian

      for ( int i = 0; i < t.size(); i++ )
        t[i] = r[i] - xm[i];
607

Francois Gygi committed
608 609 610
      x = s;
      x.herk('l','c',-0.5,t,1.0);
      x.symmetrize('l');
611

Francois Gygi committed
612 613
      for ( int i = 0; i < t.size(); i++ )
        t[i] = x[i] - xm[i];
614

Francois Gygi committed
615
      diff = t.nrm2();
616

Francois Gygi committed
617 618 619 620 621 622
      xm = x;
      iter++;
    }
    c_.hemm('r','l',1.0,x,sd.c(),1.0);
  }
}
623

624 625 626
////////////////////////////////////////////////////////////////////////////////
void SlaterDet::lowdin(void)
{
627
  cleanup();
628
  if ( basis_->real() )
629 630 631 632 633 634 635
  {
    ComplexMatrix c_tmp(c_);
    DoubleMatrix c_proxy(c_);
    DoubleMatrix c_tmp_proxy(c_tmp);
    DoubleMatrix l(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
636

637
    l.clear();
638
    l.syrk('l','t',2.0,c_proxy,0.0);
639
    l.syr('l',-1.0,c_proxy,0,'r');
640

641
    //cout << "SlaterDet::lowdin: A=\n" << l << endl;
642

643 644 645 646 647
    // Cholesky decomposition of A=Y^T Y
    l.potrf('l');
    // The lower triangle of l now contains the Cholesky factor of Y^T Y

    //cout << "SlaterDet::lowdin: L=\n" << l << endl;
648

649 650 651 652
    // Compute the polar decomposition of R = L^T

    x.transpose(1.0,l,0.0);
    // x now contains R
653

654
    const double tol = 1.e-6;
Francois Gygi committed
655
    const int maxiter = 3;
656
    x.polar(tol,maxiter);
657 658

    // x now contains the orthogonal polar factor U of the
659
    // polar decomposition R = UH
660

661
    //cout << " SlaterDet::lowdin: orthogonal polar factor=\n" << x << endl;
662

663 664 665
    // Compute L^-1
    l.trtri('l','n');
    // l now contains L^-1
666

667 668
    // Form the product L^-T U
    t.gemm('t','n',1.0,l,x,0.0);
669

670 671 672 673
    // Multiply c by L^-T U
    c_proxy.gemm('n','n',1.0,c_tmp_proxy,t,0.0);
  }
  else
674
  {
675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714
    // complex case
    ComplexMatrix c_tmp(c_);
    ComplexMatrix l(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());

    l.clear();
    l.herk('l','c',1.0,c_,0.0);

    //cout << "SlaterDet::lowdin: A=\n" << l << endl;

    // Cholesky decomposition of A=Y^H Y
    l.potrf('l');
    // The lower triangle of l now contains the Cholesky factor of Y^T Y

    //cout << "SlaterDet::lowdin: L=\n" << l << endl;

    // Compute the polar decomposition of R = L^T

    x.transpose(1.0,l,0.0);
    // x now contains R

    const double tol = 1.e-6;
    const int maxiter = 3;
    x.polar(tol,maxiter);

    // x now contains the unitary polar factor U of the
    // polar decomposition R = UH

    //cout << " SlaterDet::lowdin: unitary polar factor=\n" << x << endl;

    // Compute L^-1
    l.trtri('l','n');
    // l now contains L^-1

    // Form the product L^-T U
    t.gemm('c','n',1.0,l,x,0.0);

    // Multiply c by L^-T U
    c_.gemm('n','n',1.0,c_tmp,t,0.0);
715 716 717 718 719 720 721
  }
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::ortho_align(const SlaterDet& sd)
{
  // Orthogonalize *this and align with sd
722
  cleanup();
723
  if ( basis_->real() )
724 725 726 727 728 729 730
  {
    ComplexMatrix c_tmp(c_);
    DoubleMatrix c_proxy(c_);
    DoubleMatrix sdc_proxy(sd.c());
    DoubleMatrix c_tmp_proxy(c_tmp);
    DoubleMatrix l(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
731

732 733 734 735
#if TIMING
    tmap["syrk"].reset();
    tmap["syrk"].start();
#endif
736
    l.clear();
737
    l.syrk('l','t',2.0,c_proxy,0.0);
738
    l.syr('l',-1.0,c_proxy,0,'r');
739 740 741
#if TIMING
    tmap["syrk"].stop();
#endif
742

743
    //cout << "SlaterDet::ortho_align: A=\n" << l << endl;
744

745
    // Cholesky decomposition of A=Y^T Y
746 747 748 749
#if TIMING
    tmap["potrf"].reset();
    tmap["potrf"].start();
#endif
750
    l.potrf('l');
751 752 753
#if TIMING
    tmap["potrf"].stop();
#endif
754 755 756
    // The lower triangle of l now contains the Cholesky factor of Y^T Y

    //cout << "SlaterDet::ortho_align: L=\n" << l << endl;
757

758 759
    // Compute the polar decomposition of L^-1 B
    // where B = C^T sd.C
760

761
    // Compute B: store result in x
762 763 764 765
#if TIMING
    tmap["gemm"].reset();
    tmap["gemm"].start();
#endif
766 767 768 769
    // factor -2.0 in next line: G and -G
    x.gemm('t','n',2.0,c_proxy,sdc_proxy,0.0);
    // rank-1 update using first row of sdc_proxy() and c_proxy
    x.ger(-1.0,c_proxy,0,sdc_proxy,0);
770 771 772
#if TIMING
    tmap["gemm"].stop();
#endif
773

774 775 776
    // Form the product L^-1 B, store result in x
    // triangular solve: L X = B
    // trtrs: solve op(*this) * X = Z, output in Z
777 778 779 780
#if TIMING
    tmap["trtrs"].reset();
    tmap["trtrs"].start();
#endif
781
    l.trtrs('l','n','n',x);
782 783 784
#if TIMING
    tmap["trtrs"].stop();
#endif
785 786
    // x now contains L^-1 B

787 788 789 790 791
    // compute the polar decomposition of X = L^-1 B
#if TIMING
    tmap["polar"].reset();
    tmap["polar"].start();
#endif
792 793 794
    const double tol = 1.e-6;
    const int maxiter = 2;
    x.polar(tol,maxiter);
795 796 797
#if TIMING
    tmap["polar"].stop();
#endif
798 799

    // x now contains the orthogonal polar factor X of the
800
    // polar decomposition L^-1 B = XH
801 802

    //cout << " SlaterDet::ortho_align: orthogonal polar factor=\n"
803
    //     << x << endl;
804

805 806
    // Form the product L^-T Q
    // Solve trans(L) Z = X
807 808 809 810
#if TIMING
    tmap["trtrs2"].reset();
    tmap["trtrs2"].start();
#endif
811
    l.trtrs('l','t','n',x);
812 813 814
#if TIMING
    tmap["trtrs2"].stop();
#endif
815

816
    // x now contains L^-T Q
817

818
    // Multiply c by L^-T Q
819 820 821 822
#if TIMING
    tmap["gemm2"].reset();
    tmap["gemm2"].start();
#endif
823
    c_proxy.gemm('n','n',1.0,c_tmp_proxy,x,0.0);
824 825 826
#if TIMING
    tmap["gemm2"].stop();
#endif
827 828
  }
  else
829
  {
830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925
    // complex case
    ComplexMatrix c_tmp(c_);
    const ComplexMatrix& sdc = sd.c();
    ComplexMatrix l(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());

#if TIMING
    tmap["herk"].reset();
    tmap["herk"].start();
#endif
    l.clear();
    l.herk('l','c',1.0,c_,0.0);
#if TIMING
    tmap["herk"].stop();
#endif

    // Cholesky decomposition of A=Y^H Y
#if TIMING
    tmap["potrf"].reset();
    tmap["potrf"].start();
#endif
    l.potrf('l');
#if TIMING
    tmap["potrf"].stop();
#endif
    // The lower triangle of l now contains the Cholesky factor of Y^T Y

    //cout << "SlaterDet::ortho_align: L=\n" << l << endl;

    // Compute the polar decomposition of L^-1 B
    // where B = C^H sd.C

    // Compute B: store result in x
#if TIMING
    tmap["gemm"].reset();
    tmap["gemm"].start();
#endif
    x.gemm('c','n',1.0,c_,sdc,0.0);
#if TIMING
    tmap["gemm"].stop();
#endif

    // Form the product L^-1 B, store result in x
    // triangular solve: L X = B
    // trtrs: solve op(*this) * X = Z, output in Z
#if TIMING
    tmap["trtrs"].reset();
    tmap["trtrs"].start();
#endif
    l.trtrs('l','n','n',x);
#if TIMING
    tmap["trtrs"].stop();
#endif
    // x now contains L^-1 B

    // compute the polar decomposition of X = L^-1 B
#if TIMING
    tmap["polar"].reset();
    tmap["polar"].start();
#endif
    const double tol = 1.e-6;
    const int maxiter = 2;
    x.polar(tol,maxiter);
#if TIMING
    tmap["polar"].stop();
#endif

    // x now contains the unitary polar factor X of the
    // polar decomposition L^-1 B = XH

    //cout << " SlaterDet::ortho_align: unitary polar factor=\n"
    //     << x << endl;

    // Form the product L^-T Q
    // Solve trans(L) Z = X
#if TIMING
    tmap["trtrs2"].reset();
    tmap["trtrs2"].start();
#endif
    l.trtrs('l','c','n',x);
#if TIMING
    tmap["trtrs2"].stop();
#endif

    // x now contains L^-H Q

    // Multiply c by L^-H Q
#if TIMING
    tmap["gemm2"].reset();
    tmap["gemm2"].start();
#endif
    c_.gemm('n','n',1.0,c_tmp,x,0.0);
#if TIMING
    tmap["gemm2"].stop();
#endif

926
  }
Francois Gygi committed
927
#if TIMING
928 929 930 931 932 933 934 935 936
  for ( TimerMap::iterator i = tmap.begin(); i != tmap.end(); i++ )
  {
    double time = (*i).second.real();
    double tmin = time;
    double tmax = time;
    ctxt_.dmin(1,1,&tmin,1);
    ctxt_.dmax(1,1,&tmax,1);
    if ( ctxt_.onpe0() )
    {
937
      cout << "<timing name=\""
938 939 940
           << (*i).first << "\""
           << " min=\"" << setprecision(3) << tmin << "\""
           << " max=\"" << setprecision(3) << tmax << "\"/>"
941
           << endl;
942 943
    }
  }
Francois Gygi committed
944
#endif
945 946 947 948 949 950
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::align(const SlaterDet& sd)
{
  // Align *this with sd
951
  if ( basis_->real() )
952 953 954 955 956
  {
    ComplexMatrix c_tmp(c_);
    DoubleMatrix c_proxy(c_);
    DoubleMatrix sdc_proxy(sd.c());
    DoubleMatrix c_tmp_proxy(c_tmp);
957

958 959
    DoubleMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
960

961
    // Compute the polar decomposition of B
962
    // where B = C^H sd.C
963

Francois Gygi committed
964 965 966
#if TIMING
    tmap["align_gemm1"].start();
#endif
967 968 969 970 971
    // Compute B: store result in x
    // factor -2.0 in next line: G and -G
    x.gemm('t','n',2.0,c_proxy,sdc_proxy,0.0);
    // rank-1 update using first row of sdc_proxy() and c_proxy
    x.ger(-1.0,c_proxy,0,sdc_proxy,0);
Francois Gygi committed
972 973 974
#if TIMING
    tmap["align_gemm1"].stop();
#endif
975

976 977 978
    // x now contains B

    //cout << "SlaterDet::align: B=\n" << x << endl;
979

980
    // Compute the distance | c - sdc | before alignment
981 982
    //for ( int i = 0; i < c_proxy.size(); i++ )
    //  c_tmp_proxy[i] = c_proxy[i] - sdc_proxy[i];
Francois Gygi committed
983 984
    //cout << " SlaterDet::align: distance before: "
    //     << c_tmp_proxy.nrm2() << endl;
985

986
    // compute the polar decomposition of B
987
    double tol = 1.e-6;
Francois Gygi committed
988 989
    const int maxiter = 3;
#if TIMING
990
    tmap["align_polar"].start();
Francois Gygi committed
991
#endif
992
    x.polar(tol,maxiter);
Francois Gygi committed
993 994 995
#if TIMING
    tmap["align_while"].stop();
#endif
996 997

    // x now contains the orthogonal polar factor X of the
998
    // polar decomposition B = XH
999

1000
    //cout << " SlaterDet::align: orthogonal polar factor=\n" << x << endl;
1001

Francois Gygi committed
1002 1003 1004
#if TIMING
    tmap["align_gemm2"].start();
#endif
1005 1006 1007
    // Multiply c by X
    c_tmp_proxy = c_proxy;
    c_proxy.gemm('n','n',1.0,c_tmp_proxy,x,0.0);
Francois Gygi committed
1008 1009 1010
#if TIMING
    tmap["align_gemm2"].stop();
#endif
1011

1012
    // Compute the distance | c - sdc | after alignment
Francois Gygi committed
1013 1014 1015 1016
    //for ( int i = 0; i < c_proxy.size(); i++ )
    //  c_tmp_proxy[i] = c_proxy[i] - sdc_proxy[i];
    //cout << " SlaterDet::align: distance after:  "
    //     << c_tmp_proxy.nrm2() << endl;
1017

1018 1019
  }
  else
1020
  {
1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080
    // complex case
    ComplexMatrix c_tmp(c_);
    const ComplexMatrix& sdc = sd.c();

    ComplexMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    ComplexMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());

    // Compute the polar decomposition of B
    // where B = C^H sd.C

#if TIMING
    tmap["align_gemm1"].start();
#endif
    // Compute B: store result in x
    x.gemm('c','n',1.0,c_,sdc,0.0);
#if TIMING
    tmap["align_gemm1"].stop();
#endif

    // x now contains B

    //cout << "SlaterDet::align: B=\n" << x << endl;

    // Compute the distance | c - sdc | before alignment
    //for ( int i = 0; i < c_proxy.size(); i++ )
    //  c_tmp_proxy[i] = c_proxy[i] - sdc_proxy[i];
    //cout << " SlaterDet::align: distance before: "
    //     << c_tmp_proxy.nrm2() << endl;

    // compute the polar decomposition of B
    double tol = 1.e-6;
    const int maxiter = 3;
#if TIMING
    tmap["align_polar"].start();
#endif
    x.polar(tol,maxiter);
#if TIMING
    tmap["align_while"].stop();
#endif

    // x now contains the unitary polar factor X of the
    // polar decomposition B = XH

    //cout << " SlaterDet::align: unitary polar factor=\n" << x << endl;

#if TIMING
    tmap["align_gemm2"].start();
#endif
    // Multiply c by X
    c_tmp = c_;
    c_.gemm('n','n',1.0,c_tmp,x,0.0);
#if TIMING
    tmap["align_gemm2"].stop();
#endif

    // Compute the distance | c - sdc | after alignment
    //for ( int i = 0; i < c_proxy.size(); i++ )
    //  c_tmp_proxy[i] = c_proxy[i] - sdc_proxy[i];
    //cout << " SlaterDet::align: distance after:  "
    //     << c_tmp_proxy.nrm2() << endl;
1081 1082 1083
  }
}

1084
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
1085
complex<double> SlaterDet::dot(const SlaterDet& sd) const
1086
{
Francois Gygi committed
1087 1088
  // dot product of Slater determinants at the same kpoint: dot = tr (V^T W)
  assert(basis_->kpoint() == sd.kpoint());
1089
  if ( basis_->real() )
1090 1091 1092 1093 1094 1095
  {
    // DoubleMatrix proxy for c_ and sd.c()
    const DoubleMatrix c_proxy(c_);
    const DoubleMatrix sdc_proxy(sd.c());
    // factor 2.0: G and -G
    double d = 2.0 * c_proxy.dot(sdc_proxy);
1096

1097 1098 1099 1100 1101 1102 1103 1104 1105 1106
    // correct double counting of first element
    double sum = 0.0;
    if ( ctxt_.myrow() == 0 )
    {
      // compute the scalar product of the first rows of c_ and sd.c_
      const double *c = c_proxy.cvalptr(0);
      const double *sdc = sdc_proxy.cvalptr(0);
      int len = c_proxy.nloc();
      // stride of scalar product is mloc
      int stride = c_proxy.mloc();
Francois Gygi committed
1107
      sum = ddot(&len,c,&stride,sdc,&stride);
1108 1109 1110 1111 1112 1113
    }
    ctxt_.dsum(1,1,&sum,1);
    return d - sum;
  }
  else
  {
Francois Gygi committed
1114
    return c_.dot(sd.c());
1115 1116 1117
  }
}

Francois Gygi committed
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 1150 1151 1152
////////////////////////////////////////////////////////////////////////////////
void SlaterDet::update_occ(int nel, int nspin)
{
  // compute occupation numbers as 0.0, 1.0 or 2.0
  // if nspin = 1: use 0, 1 or 2
  // if nspin = 2: use 0 or 1;
  assert (nel >= 0);
  assert (occ_.size() == c_.n());
  if ( nspin == 1 )
  {
    assert (nel <= 2*c_.n());
    int ndouble = nel/2;
    for ( int n = 0; n < ndouble; n++ )
      occ_[n] = 2.0;
    for ( int n = ndouble; n < ndouble+nel%2; n++ )
      occ_[n] = 1.0;
    for ( int n = ndouble+nel%2; n < c_.n(); n++ )
      occ_[n] = 0.0;
  }
  else if ( nspin == 2 )
  {
    assert (nel <= c_.n());
    for ( int n = 0; n < nel; n++ )
      occ_[n] = 1.0;
    for ( int n = nel; n < c_.n(); n++ )
      occ_[n] = 0.0;
  }
  else
  {
    // incorrect value of nspin_
    assert(false);
  }
}

////////////////////////////////////////////////////////////////////////////////
1153
double SlaterDet::total_charge(void) const
Francois Gygi committed
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
{
  // compute total charge from occ_[i]
  double sum = 0.0;
  for ( int n = 0; n < occ_.size(); n++ )
  {
    sum += occ_[n];
  }
  return sum;
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::update_occ(int nspin, double mu, double temp)
{
  // compute occupation numbers using a Fermi distribution f(mu,temp)
  // and the eigenvalues in eig_[i]
  assert(nspin==1 || nspin==2);
  assert (occ_.size() == c_.n());
  assert (eig_.size() == c_.n());
  if ( nspin == 1 )
  {
    for ( int n = 0; n < eig_.size(); n++ )
    {
      occ_[n] = 2.0 * fermi(eig_[n],mu,temp);
    }
  }
  else if ( nspin == 2 )
  {
    for ( int n = 0; n < eig_.size(); n++ )
    {
      occ_[n] = fermi(eig_[n],mu,temp);
    }
  }
  else
  {
    // incorrect value of nspin_
    assert(false);
  }
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::fermi(double e, double mu, double fermitemp)
{
  // e, mu in Hartree, fermitemp in Kelvin

  if ( fermitemp == 0.0 )
  {
    if ( e < mu ) return 1.0;
    else if ( e == mu ) return 0.5;
    else return 0.0;
  }
  const double kb = 3.1667907e-6; // Hartree/Kelvin
  const double kt = kb * fermitemp;
  double arg = ( e - mu ) / kt;

  if ( arg < -30.0 ) return 1.0;
  if ( arg >  30.0 ) return 0.0;

  return 1.0 / ( 1.0 + exp ( arg ) );
}

////////////////////////////////////////////////////////////////////////////////
1215
double SlaterDet::entropy(int nspin) const
Francois Gygi committed
1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234
{
  // return dimensionless entropy
  // the contribution to the free energy is - t_kelvin * k_boltz * wf.entropy()

  assert(nspin==1 || nspin==2);
  const double fac = ( nspin > 1 ) ? 1.0 : 2.0;
  double sum = 0.0;
  for ( int n = 0; n < occ_.size(); n++ )
  {
    const double f = occ_[n] / fac;
    if ( f > 0.0  &&  f < 1.0 )
    {
      sum -= fac * ( f * log(f) + (1.0-f) * log(1.0-f) );
    }
  }
  return sum;
}

////////////////////////////////////////////////////////////////////////////////
1235
double SlaterDet::ortho_error(void) const
Francois Gygi committed
1236 1237 1238
{
  // deviation from orthogonality of c_
  double error;
1239
  if ( basis_->real() )
Francois Gygi committed
1240 1241 1242 1243
  {
    // k = 0 case
    // declare a proxy DoubleMatrix for c_
    DoubleMatrix c_proxy(c_);
1244

Francois Gygi committed
1245
    DoubleMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
1246

Francois Gygi committed
1247 1248 1249
    // real symmetric rank-k update
    // factor 2.0 in next line: G and -G
    s.syrk('l','t',2.0,c_proxy,0.0); // compute real overlap matrix
1250

Francois Gygi committed
1251 1252 1253
    // correct for double counting of G=0
    // symmetric rank-1 update using first row of c_proxy
    s.syr('l',-1.0,c_proxy,0,'r');
1254

Francois Gygi committed
1255 1256
    DoubleMatrix id(ctxt_,s.m(),s.n(),s.mb(),s.nb());
    id.identity();
1257

Francois Gygi committed
1258
    s -= id; // subtract identity matrix from S
1259

Francois Gygi committed
1260 1261 1262 1263 1264
    error = s.nrm2();
  }
  else
  {
    // k != 0 case
1265

Francois Gygi committed
1266 1267
    ComplexMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    s.herk('l','c',1.0,c_,0.0);
1268

Francois Gygi committed
1269 1270
    ComplexMatrix id(ctxt_,s.m(),s.n(),s.mb(),s.nb());
    id.identity();
1271

Francois Gygi committed
1272
    s -= id; // subtract identity matrix from S
1273

Francois Gygi committed
1274 1275 1276 1277 1278 1279 1280 1281
    error = s.nrm2();
  }
  return error;
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::randomize(double amplitude)
{
1282 1283
  if ( basis_->size() == 0 )
    return;
1284

Francois Gygi committed
1285 1286 1287
  for ( int n = 0; n < c_.nloc(); n++ )
  {
    complex<double>* p = c_.valptr(c_.mloc()*n);
1288
    for ( int i = 0; i < basis_->localsize(); i++ )
Francois Gygi committed
1289 1290 1291 1292 1293 1294
    {
      double re = drand48();
      double im = drand48();
      p[i] += amplitude * complex<double>(re,im);
    }
  }
1295
  // gram does an initial cleanup
Francois Gygi committed
1296 1297 1298 1299 1300 1301
  gram();
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::cleanup(void)
{
1302
  // set Im( c(G=0) ) to zero for real case and
Francois Gygi committed
1303
  // set the empty rows of the matrix c_ to zero
1304 1305 1306
  // The empty rows are located between i = basis_->localsize() and
  // c_.mloc(). Empty rows are necessary to insure that the
  // local size c_.mloc() is the same on all processes, while the
Francois Gygi committed
1307 1308 1309 1310 1311
  // local basis size is not.
  for ( int n = 0; n < c_.nloc(); n++ )
  {
    complex<double>* p = c_.valptr(c_.mloc()*n);
    // reset imaginary part of G=0 component to zero
1312
    if ( basis_->real() && c_.mloc() > 0 && ctxt_.myrow() == 0 )
Francois Gygi committed
1313 1314
    {
      // index of G=0 element
1315
      p[0] = complex<double> ( p[0].real(), 0.0);
Francois Gygi committed
1316 1317
    }
    // reset values of empty rows of c_ to zero
1318
    for ( int i = basis_->localsize(); i < c_.mloc(); i++ )
Francois Gygi committed
1319 1320 1321 1322 1323 1324 1325 1326 1327 1328
      p[i] = 0.0;
  }
}

////////////////////////////////////////////////////////////////////////////////
SlaterDet& SlaterDet::operator=(SlaterDet& rhs)
{
  if ( this == &rhs ) return *this;
  assert(ctxt_.ictxt() == rhs.context().ictxt());
  c_ = rhs.c_;
1329 1330
  occ_ = rhs.occ_;
  eig_ = rhs.eig_;
Francois Gygi committed
1331 1332 1333 1334 1335 1336
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::memsize(void) const
{
1337
  return basis_->memsize() + c_.memsize();
Francois Gygi committed
1338 1339 1340 1341 1342
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::localmemsize(void) const
{
1343
  return basis_->localmemsize() + c_.localmemsize();
Francois Gygi committed
1344 1345 1346
}

////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
1347
void SlaterDet::print(ostream& os, string encoding, double weight, int ispin,
1348
  int nspin) const
Francois Gygi committed
1349
{
1350
  FourierTransform ft(*basis_,basis_->np(0),basis_->np(1),basis_->np(2));
Francois Gygi committed
1351
  vector<complex<double> > wftmp(ft.np012loc());
Francois Gygi committed
1352
  const bool real_basis = basis_->real();
1353
  const int wftmpr_size = real_basis ? ft.np012() : 2*ft.np012();
1354
  const int wftmpr_loc_size = real_basis ? ft.np012loc() : 2*ft.np012loc();
Francois Gygi committed
1355
  vector<double> wftmpr(wftmpr_size);
Francois Gygi committed
1356
  Base64Transcoder xcdr;
1357

Francois Gygi committed
1358 1359
  if ( ctxt_.onpe0() )
  {
Francois Gygi committed
1360 1361 1362 1363 1364
    string spin = (ispin > 0) ? "down" : "up";
    os << "<slater_determinant";
    if ( nspin == 2 )
      os << " spin=\"" << spin << "\"";
    os << " kpoint=\"" << basis_->kpoint() << "\"\n"
1365
       << "  weight=\"" << setprecision(12) <<  weight << "\""
Francois Gygi committed
1366
       << " size=\"" << nst() << "\">" << endl;
1367 1368

    os << "<density_matrix form=\"diagonal\" size=\"" << nst() << "\">"
Francois Gygi committed
1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381
       << endl;
    os.setf(ios::fixed,ios::floatfield);
    os.setf(ios::right,ios::adjustfield);
    for ( int i = 0; i < nst(); i++ )
    {
      os << " " << setprecision(8) << occ_[i];
      if ( i%10 == 9 )
        os << endl;
    }
    if ( nst()%10 != 0 )
      os << endl;
    os << "</density_matrix>" << endl;
  }
1382

Francois Gygi committed
1383 1384
  for ( int n = 0; n < nst(); n++ )
  {
1385
    // Barrier to limit the number of messages sent to task 0
1386 1387
    // that don't have a receive posted
    ctxt_.barrier();
1388

1389
    // transform data on ctxt_.mycol()
Francois Gygi committed
1390 1391
    if ( c_.pc(n) == ctxt_.mycol() )
    {
1392
      //cout << " state " << n << " is stored on column "
Francois Gygi committed
1393 1394 1395
      //     << ctxt_.mycol() << " local index: " << c_.y(n) << endl;
      int nloc = c_.y(n); // local index
      ft.backward(c_.cvalptr(c_.mloc()*nloc),&wftmp[0]);
1396

Francois Gygi committed
1397 1398 1399 1400 1401 1402 1403 1404
      if ( real_basis )
      {
        double *a = (double*) &wftmp[0];
        for ( int i = 0; i < ft.np012loc(); i++ )
          wftmpr[i] = a[2*i];
      }
      else
      {
1405 1406
        memcpy((void*)&wftmpr[0],(void*)&wftmp[0],
               ft.np012loc()*sizeof(complex<double>));
Francois Gygi committed
1407
      }
1408
    }
1409

1410 1411 1412 1413
    // send blocks of wftmpr to pe0
    for ( int i = 0; i < ctxt_.nprow(); i++ )
    {
      bool iamsending = c_.pc(n) == ctxt_.mycol() && i == ctxt_.myrow();
1414

1415 1416 1417
      // send size of wftmpr block
      int size=-1;
      if ( ctxt_.onpe0() )
Francois Gygi committed
1418
      {
1419
        if ( iamsending )
Francois Gygi committed
1420
        {
1421 1422 1423 1424 1425 1426 1427 1428 1429
          // sending to self, size not needed
        }
        else
          ctxt_.irecv(1,1,&size,1,i,c_.pc(n));
      }
      else
      {
        if ( iamsending )
        {
1430
          size = wftmpr_loc_size;
Francois Gygi committed
1431
          ctxt_.isend(1,1,&size,1,0,0);
1432 1433
        }
      }
1434

1435 1436 1437 1438 1439 1440 1441 1442 1443 1444
      // send wftmpr block
      if ( ctxt_.onpe0() )
      {
        if ( iamsending )
        {
          // do nothing, data is already in place
        }
        else
        {
          int istart = ft.np0() * ft.np1() * ft.np2_first(i);
Francois Gygi committed
1445 1446
          if ( !real_basis )
            istart *= 2;
1447
          ctxt_.drecv(size,1,&wftmpr[istart],size,i,c_.pc(n));
1448 1449 1450 1451 1452 1453
        }
      }
      else
      {
        if ( iamsending )
        {
1454
          ctxt_.dsend(size,1,&wftmpr[0],size,0,0);
Francois Gygi committed
1455 1456 1457
        }
      }
    }
1458

1459
    // process the data
Francois Gygi committed
1460 1461 1462
    if ( ctxt_.onpe0() )
    {
      // wftmpr is now complete on task 0
Francois Gygi committed
1463 1464 1465
      // wftmpr contains either a real of a complex array

      const string element_t