SlaterDet.C 39.8 KB
Newer Older
Francois Gygi committed
1 2 3 4 5
////////////////////////////////////////////////////////////////////////////////
//
// SlaterDet.C
//
////////////////////////////////////////////////////////////////////////////////
6
// $Id: SlaterDet.C,v 1.38 2006-05-29 01:15:36 fgygi Exp $
Francois Gygi committed
7 8 9 10 11

#include "SlaterDet.h"
#include "FourierTransform.h"
#include "Context.h"
#include "blas.h" // daxpy
Francois Gygi committed
12
#include "Base64Transcoder.h"
Francois Gygi committed
13
#include "Timer.h"
Francois Gygi committed
14 15 16 17 18 19 20 21 22 23 24 25

#include <cstdlib>
#include <iostream>
#include <iomanip>
#if USE_CSTDIO_LFS
#include <sstream>
#include <cstdio>
#endif
using namespace std;

////////////////////////////////////////////////////////////////////////////////
SlaterDet::SlaterDet(const Context& ctxt, D3vector kpoint) : ctxt_(ctxt), 
26
 c_(ctxt)
Francois Gygi committed
27 28
{
  //cout << ctxt.mype() << ": SlaterDet::SlaterDet: ctxt.mycol="
29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
  //     << ctxt.mycol() << " basis_->context(): " 
  //     << basis_->context();
  my_col_ctxt_ = 0;
  for ( int icol = 0; icol < ctxt_.npcol(); icol++ )
  {
    Context* col_ctxt = new Context(ctxt_,ctxt_.nprow(),1,0,icol);
    ctxt_.barrier();
    if ( icol == ctxt_.mycol() )
      my_col_ctxt_ = col_ctxt;
    else
      delete col_ctxt;
  }
  //cout << ctxt_.mype() << ": SlaterDet::SlaterDet: my_col_ctxt: " 
  //     << *my_col_ctxt_;
  basis_ = new Basis(*my_col_ctxt_,kpoint);
Francois Gygi committed
44
}
Francois Gygi committed
45 46 47

////////////////////////////////////////////////////////////////////////////////
SlaterDet::SlaterDet(const SlaterDet& rhs) : ctxt_(rhs.context()),
48 49
  basis_(new Basis(*(rhs.basis_))), 
  my_col_ctxt_(new Context(*(rhs.my_col_ctxt_))) , c_(rhs.c_){}
Francois Gygi committed
50 51 52
  
////////////////////////////////////////////////////////////////////////////////
SlaterDet::~SlaterDet() 
53 54 55
{
  delete basis_;
  delete my_col_ctxt_;
Francois Gygi committed
56 57 58 59 60 61 62
  // cout << ctxt_.mype() << ": SlaterDet::dtor: ctxt=" << ctxt_;
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::resize(const UnitCell& cell, const UnitCell& refcell, 
  double ecut, int nst)
{
63 64
  //!! Test in next line should be replaced by test on basis min/max indices
  //!! to signal change in basis vectors
65
  //if ( basis_->refcell().volume() != 0.0 && !refcell.encloses(cell) )
66 67 68
  //{
    //!! << " SlaterDet::resize: warning" << endl;
    //cout << " SlaterDet::resize: cell=" << cell;
69
    //cout << " SlaterDet::resize: refcell=" << basis_->refcell();
70
    //throw SlaterDetException("could not resize: cell not in refcell");
71
  //}
72
  
Francois Gygi committed
73 74
  try
  {
75
    basis_->resize(cell,refcell,ecut);
Francois Gygi committed
76 77 78
    occ_.resize(nst);
    eig_.resize(nst);
    
79
    const int mb = basis_->maxlocalsize();
80 81 82 83 84 85 86 87 88 89 90 91
    const int m = ctxt_.nprow() * mb;
    const int nb = nst/ctxt_.npcol() + (nst%ctxt_.npcol() > 0 ? 1 : 0);
    
    // Determine if plane wave coefficients must be reset after the resize
    // This is needed if the dimensions of the matrix c_ must be changed
    const bool needs_reset = 
      m!=c_.m() || nst!=c_.n() || mb!=c_.mb() || nb!=c_.nb();
  
    c_.resize(m,nst,mb,nb);
    
    if ( needs_reset )
      reset();
Francois Gygi committed
92 93 94 95 96 97 98
  }
  catch ( bad_alloc )
  {
    cout << " bad_alloc exception caught in SlaterDet::resize" << endl;
    throw;
  }
}
99 100 101 102
////////////////////////////////////////////////////////////////////////////////
void SlaterDet::reset(void)
{
  // initialize coefficients with lowest plane waves
103
  if ( c_.n() <= basis_->size() )
104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119
  {
    // initialize c_
    c_.clear();
    const double s2i = 1.0 / sqrt(2.0);
 
    // 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;
 
      // find process row holding the smallest g vector
120
      double g2 = basis_->g2(basis_->isort(ismallest));
121 122 123 124 125 126 127 128 129 130
      // cout << "smallest vector on proc " << ctxt_.mype()
      //      << " has norm " << g2 << endl;
      int minrow, mincol;
      ctxt_.dmin('c',' ',1,1,&g2,1,&minrow,&mincol,1,-1,-1);
 
      // find column hosting state n
      int pc = c_.pc(n);
      int pr = minrow;
      if ( pr == ctxt_.myrow() )
      {
131
        int iii = basis_->isort(ismallest);
132 133 134 135 136
        ismallest++; // increment on entire process row
        if ( pc == ctxt_.mycol() )
        {
          // cout << " n=" << n << " on process "
          //      << pr << "," << pc
137 138 139 140
          //      << " vector " << basis_->idx(3*iii) << " "
          //      << basis_->idx(3*iii+1) << " "
          //      << basis_->idx(3*iii+2) << " norm="
          //      << basis_->g2(iii) << " "
141 142 143 144 145 146 147 148 149
          //      << 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
150 151 152 153 154

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::compute_density(FourierTransform& ft, 
  double weight, double* rho) const
{
Francois Gygi committed
155
  //Timer tm_ft, tm_rhosum;
Francois Gygi committed
156 157 158 159
  // compute density of the states residing on my column of ctxt_
  assert(occ_.size() == c_.n());
  vector<complex<double> > tmp(ft.np012loc());
  
160 161
  assert(basis_->cell().volume() > 0.0);
  const double omega_inv = 1.0 / basis_->cell().volume();
Francois Gygi committed
162 163 164 165 166
  const int np012loc = ft.np012loc();
  
  for ( int i = 0; i < np012loc; i++ )
    rho[i] = 0.0;

167
  if ( basis_->real() )
Francois Gygi committed
168 169 170 171 172 173 174 175 176 177 178
  {
    // 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;
      const double fac1 = omega_inv * occ_[nn];
      const double fac2 = omega_inv * occ_[nn+1];
      
      if ( fac1 + fac2 > 0.0 )
      {
Francois Gygi committed
179
        //tm_ft.start();
Francois Gygi committed
180
        ft.backward(c_.cvalptr(n*c_.mloc()),
Francois Gygi committed
181 182
                    c_.cvalptr((n+1)*c_.mloc()),&tmp[0]);
        //tm_ft.stop();
Francois Gygi committed
183 184
        const double* psi = (double*) &tmp[0];
        int ii = 0;
Francois Gygi committed
185
        //tm_rhosum.start();
Francois Gygi committed
186 187 188 189 190 191 192
        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
193
        //tm_rhosum.start();
Francois Gygi committed
194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233
      }
    }
    if ( nstloc() % 2 != 0 )
    {
      const int n = nstloc()-1;
      // global n index
      const int nn = ctxt_.mycol() * c_.nb() + n;
      const double fac1 = omega_inv * occ_[nn];
      
      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;
      const double fac = omega_inv * occ_[nn];
      
      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]);
      }
    }
  }
Francois Gygi committed
234 235 236 237
  // cout << "SlaterDet: compute_density: ft_bwd time: " 
  //      << tm_ft.real() << endl;
  // cout << "SlaterDet: compute_density: rhosum time: " 
  //      << tm_rhosum.real() << endl;
Francois Gygi committed
238 239 240 241
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::rs_mul_add(FourierTransform& ft, 
Francois Gygi committed
242
  const double* v, SlaterDet& sdp) const
Francois Gygi committed
243 244 245 246 247 248 249 250 251 252 253 254 255
{
  // 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]
  
  vector<complex<double> > tmp(ft.np012loc());
  vector<complex<double> > ctmp(2*c_.mloc());
  
  const int np012loc = ft.np012loc();
  const int mloc = c_.mloc();
  double* p = (double*) &tmp[0];
  double* dcp = (double*) sdp.c().valptr();

256
  if ( basis_->real() )
Francois Gygi committed
257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276
  {
    // 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]);
      int ii = 0;
      for ( int i = 0; i < np012loc; i++ )
      {
        const double psi1 = p[ii];
        const double psi2 = p[ii+1];
        const double vii = v[i];
        p[ii]   = vii * psi1;
        p[ii+1] = vii * psi2;
        ii++; ii++;
      }
      ft.forward(&tmp[0], &ctmp[0], &ctmp[mloc]);
      int len = 4 * mloc;
      int inc1 = 1;
      double alpha = 1.0;
Francois Gygi committed
277
      daxpy(&len,&alpha,(double*)&ctmp[0],&inc1,&dcp[2*n*mloc],&inc1);
Francois Gygi committed
278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295
    }
    if ( nstloc() % 2 != 0 )
    {
      const int n = nstloc()-1;
      ft.backward(c_.cvalptr(n*mloc),&tmp[0]);
      int ii = 0;
      for ( int i = 0; i < np012loc; i++ )
      {
        const double psi1 = p[ii];
        const double vii = v[i];
        p[ii]   = vii * psi1;
        p[ii+1] = 0.0;
        ii++; ii++;
      }
      ft.forward(&tmp[0], &ctmp[0]);
      int len = 2 * mloc;
      int inc1 = 1;
      double alpha = 1.0;
Francois Gygi committed
296
      daxpy(&len,&alpha,(double*)&ctmp[0],&inc1,&dcp[2*n*mloc],&inc1);
Francois Gygi committed
297 298 299 300 301 302 303 304 305 306 307 308 309 310
    }
  }
  else
  {
    // only one transform at a time
    for ( int n = 0; n < nstloc(); n++ )
    {
      ft.backward(c_.cvalptr(n*mloc),&tmp[0]);
      for ( int i = 0; i < np012loc; i++ )
        tmp[i] *= v[i];
      ft.forward(&tmp[0], &ctmp[0]);
      int len = 2 * mloc;
      int inc1 = 1;
      double alpha = 1.0;
Francois Gygi committed
311
      daxpy(&len,&alpha,(double*)&ctmp[0],&inc1,&dcp[2*n*mloc],&inc1);
Francois Gygi committed
312 313 314 315 316 317 318 319
    }
  }
  
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::gram(void)
{
320
  if ( basis_->real() )
Francois Gygi committed
321 322 323 324 325 326 327
  {
    // k = 0 case
    // create a DoubleMatrix proxy for c_
    DoubleMatrix c_proxy(c_);
    DoubleMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    s.syrk('l','t',2.0,c_proxy,0.0); 
    s.syr('l',-1.0,c_proxy,0,'r');
328 329
#ifdef CHOLESKY_REMAP
    // create a square context for the Cholesky decomposition
330 331
    // int nsq = (int) sqrt((double) ctxt_.size());
    int nsq = CHOLESKY_REMAP;
332 333 334 335 336
    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);
337 338
#else
    s.potrf('l'); // Cholesky decomposition: S = L * L^T
339
#endif
Francois Gygi committed
340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356
    // solve triangular system X * L^T = C
    c_proxy.trsm('r','l','t','n',1.0,s);
  }
  else
  {
    // k != 0 case
    ComplexMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    s.herk('l','c',1.0,c_,1.0);
    s.potrf('l'); // Cholesky decomposition: S = L * L^T
    // solve triangular system X * L^T = C
    c_.trsm('r','l','c','n',1.0,s);
  }
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::riccati(SlaterDet& sd)
{
357
  if ( basis_->real() )
Francois Gygi committed
358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 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 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469
  {
    // 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();
 
    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());
    
    // DoubleMatrix proxy for c_ and sd.c()
    DoubleMatrix c_proxy(c_);
    DoubleMatrix sdc_proxy(sd.c());
    
    // 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');
    // 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);
 
    xm = s;
    xm.symmetrize('l');
 
    s.syrk('l','t',0.5,r,1.0); // s = s + 0.5 * r^T * r
    s.symmetrize('l');
 
    double diff = 1.0;
    const double epsilon = 1.e-10;
    const int maxiter = 20;
    int iter = 0;
 
    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];
 
      x = s;
      x.syrk('l','t',-0.5,t,1.0);
 
      // get full matrix x
      x.symmetrize('l');
 
      for ( int i = 0; i < t.size(); i++ )
        t[i] = x[i] - xm[i];
 
      diff = t.nrm2();
 
      xm = x;
      iter++;
    }
    c_proxy.symm('r','l',1.0,x,sdc_proxy,1.0);
  }
  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();
 
    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());
    
    // 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);
 
    xm = s;
    xm.symmetrize('l');
 
    // s = s + 0.5 * r^H * r
    s.herk('l','c',0.5,r,1.0);
    s.symmetrize('l');
 
    double diff = 1.0;
    const double epsilon = 1.e-10;
    const int maxiter = 20;
    int iter = 0;
 
    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];
 
      x = s;
      x.herk('l','c',-0.5,t,1.0);
      x.symmetrize('l');
 
      for ( int i = 0; i < t.size(); i++ )
        t[i] = x[i] - xm[i];
 
      diff = t.nrm2();
 
      xm = x;
      iter++;
    }
    c_.hemm('r','l',1.0,x,sd.c(),1.0);
  }
}
  
470 471 472 473
////////////////////////////////////////////////////////////////////////////////
void SlaterDet::lowdin(void)
{
  // Higham algorithm for polar decomposition
474
  if ( basis_->real() )
475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 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
  {
    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 xp(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    
    l.clear();
    l.syrk('l','t',2.0,c_proxy,0.0); 
    l.syr('l',-1.0,c_proxy,0,'r');
    
    //cout << "SlaterDet::lowdin: A=\n" << l << endl;
    
    // 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;
    
    // Compute the polar decomposition of R = L^T

    x.transpose(1.0,l,0.0);
    // x now contains R
    //cout << "SlaterDet::lowdin: R=\n" << x << endl;
    
    double diff = 1.0;
    const double epsilon = 1.e-10;
    const int maxiter = 20;
    int iter = 0;
 
    while ( iter < maxiter && diff > epsilon )
    {
      // t = X^T
      t.transpose(1.0,x,0.0);
      t.inverse();
      
      // t now contains X^-T
      
      // xp = 0.5 * ( x + x^-T );
      for ( int i = 0; i < x.size(); i++ )
        xp[i] = 0.5 * ( x[i] + t[i] );
      
 
      // Next lines: use t as temporary to compute || x - xp ||_F
      for ( int i = 0; i < t.size(); i++ )
        t[i] = x[i] - xp[i];
 
      diff = t.nrm2();
      
      //cout << " SlaterDet::lowdin: diff=" << diff << endl;
      
      x = xp;
      //cout << "SlaterDet::lowdin: X=\n" << x << endl;
 
      iter++;
    }
    
    // x now contains the orthogonal polar factor U of the 
    // polar decomposition R = UH
    
    //cout << " SlaterDet::lowdin: orthogonal 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('t','n',1.0,l,x,0.0);
    
    // Multiply c by L^-T U
    c_proxy.gemm('n','n',1.0,c_tmp_proxy,t,0.0);
    
  }
  else
  { 
    // complex case: not implemented
    assert(false);
  }
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::ortho_align(const SlaterDet& sd)
{
  // Orthogonalize *this and align with sd
  // Higham algorithm for polar decomposition
562
  if ( basis_->real() )
563 564 565 566 567 568 569 570 571 572
  {
    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());
    DoubleMatrix xp(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    
573 574 575 576
#if TIMING
    tmap["syrk"].reset();
    tmap["syrk"].start();
#endif
577 578 579
    l.clear();
    l.syrk('l','t',2.0,c_proxy,0.0); 
    l.syr('l',-1.0,c_proxy,0,'r');
580 581 582
#if TIMING
    tmap["syrk"].stop();
#endif
583 584 585 586
    
    //cout << "SlaterDet::ortho_align: A=\n" << l << endl;
    
    // Cholesky decomposition of A=Y^T Y
587 588 589 590
#if TIMING
    tmap["potrf"].reset();
    tmap["potrf"].start();
#endif
591
    l.potrf('l');
592 593 594
#if TIMING
    tmap["potrf"].stop();
#endif
595 596 597 598 599 600 601 602
    // 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^T sd.C
    
    // Compute B: store result in x
603 604 605 606
#if TIMING
    tmap["gemm"].reset();
    tmap["gemm"].start();
#endif
607 608 609 610
    // 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);
611 612 613
#if TIMING
    tmap["gemm"].stop();
#endif
614 615 616 617
    
    // Form the product L^-1 B, store result in x
    // triangular solve: L X = B
    // trtrs: solve op(*this) * X = Z, output in Z
618 619 620 621
#if TIMING
    tmap["trtrs"].reset();
    tmap["trtrs"].start();
#endif
622
    l.trtrs('l','n','n',x);
623 624 625
#if TIMING
    tmap["trtrs"].stop();
#endif
626 627 628 629 630 631 632 633 634 635
    // x now contains L^-1 B

    //cout << "SlaterDet::ortho_align: L^-1 B=\n" << x << endl;
    
    // compute the polar decomposition of L^-1 B
    double diff = 1.0;
    const double epsilon = 1.e-10;
    const int maxiter = 20;
    int iter = 0;
 
636 637 638 639
#if TIMING
    tmap["transpose"].reset();
    tmap["inverse"].reset();
#endif
640 641 642
    while ( iter < maxiter && diff > epsilon )
    {
      // t = X^T
643 644 645
#if TIMING
      tmap["transpose"].start();
#endif
646
      t.transpose(1.0,x,0.0);
647 648 649 650
#if TIMING
      tmap["transpose"].stop();
      tmap["inverse"].start();
#endif
651
      t.inverse();
652 653 654
#if TIMING
      tmap["inverse"].stop();
#endif
655 656 657 658 659 660 661 662 663 664 665 666
      
      // t now contains X^-T
      
      // xp = 0.5 * ( x + x^-T );
      for ( int i = 0; i < x.size(); i++ )
        xp[i] = 0.5 * ( x[i] + t[i] );
      
 
      // Next lines: use t as temporary to compute || x - xp ||_F
      for ( int i = 0; i < t.size(); i++ )
        t[i] = x[i] - xp[i];
 
667 668 669
#if TIMING
      tmap["nrm2"].start();
#endif
670
      diff = t.nrm2();
671 672 673
#if TIMING
      tmap["nrm2"].stop();
#endif
674 675 676 677 678 679 680 681 682 683 684 685
      
      //cout << " SlaterDet::ortho_align: diff=" << diff << endl;
      
      x = xp;
      //cout << "SlaterDet::ortho_align: X=\n" << x << endl;
 
      iter++;
    }
    
    // x now contains the orthogonal polar factor X of the 
    // polar decomposition L^-1 B = XH
    
686 687
    //cout << " SlaterDet::ortho_align: orthogonal polar factor=\n" 
    //     << x << endl;
688 689 690
    
    // Form the product L^-T Q
    // Solve trans(L) Z = X
691 692 693 694
#if TIMING
    tmap["trtrs2"].reset();
    tmap["trtrs2"].start();
#endif
695
    l.trtrs('l','t','n',x);
696 697 698
#if TIMING
    tmap["trtrs2"].stop();
#endif
699 700 701 702
    
    // x now contains L^-T Q
    
    // Multiply c by L^-T Q
703 704 705 706
#if TIMING
    tmap["gemm2"].reset();
    tmap["gemm2"].start();
#endif
707
    c_proxy.gemm('n','n',1.0,c_tmp_proxy,x,0.0);
708 709 710
#if TIMING
    tmap["gemm2"].stop();
#endif
711 712 713 714 715 716 717
    
  }
  else
  { 
    // complex case: not implemented
    assert(false);
  }
718 719 720 721 722 723 724 725 726 727 728 729 730 731 732
  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() )
    {
      cout << "<!-- timing "
           << setw(15) << (*i).first
           << " : " << setprecision(3) << setw(9) << tmin
           << " "   << setprecision(3) << setw(9) << tmax << " -->" << endl;
    }
  }
733 734 735 736 737 738 739
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::align(const SlaterDet& sd)
{
  // Align *this with sd
  // Higham algorithm for polar decomposition
740
  if ( basis_->real() )
741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764
  {
    ComplexMatrix c_tmp(c_);
    DoubleMatrix c_proxy(c_);
    DoubleMatrix sdc_proxy(sd.c());
    DoubleMatrix c_tmp_proxy(c_tmp);
    DoubleMatrix x(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix xp(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    DoubleMatrix t(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    
    
    // Compute the polar decomposition of B
    // where B = C^T sd.C
    
    // 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);
    
    // x now contains B

    //cout << "SlaterDet::align: B=\n" << x << endl;
    
    // Compute the distance | c - sdc | before alignment
765 766
    //for ( int i = 0; i < c_proxy.size(); i++ )
    //  c_tmp_proxy[i] = c_proxy[i] - sdc_proxy[i];
Francois Gygi committed
767 768
    //cout << " SlaterDet::align: distance before: "
    //     << c_tmp_proxy.nrm2() << endl;
769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789
    
    // compute the polar decomposition of B
    double diff = 1.0;
    const double epsilon = 1.e-10;
    const int maxiter = 20;
    int iter = 0;
 
    while ( iter < maxiter && diff > epsilon )
    {
      // t = X^T
      t.transpose(1.0,x,0.0);
      t.inverse();
      
      // t now contains X^-T
      
      // xp = 0.5 * ( x + x^-T );
      for ( int i = 0; i < x.size(); i++ )
        xp[i] = 0.5 * ( x[i] + t[i] );
      
 
      // Next lines: use t as temporary to compute || x - xp ||_F
Francois Gygi committed
790 791
      //for ( int i = 0; i < t.size(); i++ )
      //  t[i] = x[i] - xp[i];
792
 
Francois Gygi committed
793
      //diff = t.nrm2();
794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812
      
      //cout << " SlaterDet::align: diff=" << diff << endl;
      
      x = xp;
      //cout << "SlaterDet::align: X=\n" << x << endl;
 
      iter++;
    }
    
    // x now contains the orthogonal polar factor X of the 
    // polar decomposition B = XH
    
    //cout << " SlaterDet::align: orthogonal polar factor=\n" << x << endl;
        
    // Multiply c by X
    c_tmp_proxy = c_proxy;
    c_proxy.gemm('n','n',1.0,c_tmp_proxy,x,0.0);
    
    // Compute the distance | c - sdc | after alignment
Francois Gygi committed
813 814 815 816
    //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;
817 818 819 820 821 822 823 824 825
    
  }
  else
  { 
    // complex case: not implemented
    assert(false);
  }
}

826 827 828 829
////////////////////////////////////////////////////////////////////////////////
double SlaterDet::dot(const SlaterDet& sd) const
{
  // dot product of Slater determinants: dot = tr (V^T W)
830
  if ( basis_->real() )
831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847
  {
    // 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);
    
    // 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
848
      sum = ddot(&len,c,&stride,sdc,&stride);
849 850 851 852 853 854 855 856 857 858 859
    }
    ctxt_.dsum(1,1,&sum,1);
    return d - sum;
  }
  else
  {
    assert(false);
    return 0.0;
  }
}

Francois Gygi committed
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 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 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 967 968 969 970 971 972 973 974 975 976 977 978 979 980
////////////////////////////////////////////////////////////////////////////////
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);
  }
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::total_charge(void)
{
  // 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 ) );
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::entropy(int nspin)
{
  // 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;
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::ortho_error(void)
{
  // deviation from orthogonality of c_
  double error;
981
  if ( basis_->real() )
Francois Gygi committed
982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023
  {
    // k = 0 case
    // declare a proxy DoubleMatrix for c_
    DoubleMatrix c_proxy(c_);
    
    DoubleMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    
    // 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
 
    // 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');
 
    DoubleMatrix id(ctxt_,s.m(),s.n(),s.mb(),s.nb());
    id.identity();
    
    s -= id; // subtract identity matrix from S
    
    error = s.nrm2();
  }
  else
  {
    // k != 0 case
    
    ComplexMatrix s(ctxt_,c_.n(),c_.n(),c_.nb(),c_.nb());
    s.herk('l','c',1.0,c_,0.0);
    
    ComplexMatrix id(ctxt_,s.m(),s.n(),s.mb(),s.nb());
    id.identity();
 
    s -= id; // subtract identity matrix from S
    
    error = s.nrm2();
  }
  return error;
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::randomize(double amplitude)
{
1024 1025
  if ( basis_->size() == 0 )
    return;
Francois Gygi committed
1026 1027 1028 1029 1030
  // Note: randomization results depend on the process grid size and shape
  srand48(ctxt_.myproc());
  for ( int n = 0; n < c_.nloc(); n++ )
  {
    complex<double>* p = c_.valptr(c_.mloc()*n);
1031
    for ( int i = 0; i < basis_->localsize(); i++ )
Francois Gygi committed
1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046
    {
      double re = drand48();
      double im = drand48();
      p[i] += amplitude * complex<double>(re,im);
    }
  }
  cleanup();
  gram();
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::cleanup(void)
{
  // set Im( c(G=0) ) to zero and 
  // set the empty rows of the matrix c_ to zero
1047
  // The empty rows are located between i = basis_->localsize() and 
Francois Gygi committed
1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058
  // c_.mloc(). Empty rows are necessary to insure that the 
  // local size c_.mloc() is the same on all processes, while the 
  // 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
    if ( ctxt_.myrow() == 0 )
    {
      // index of G=0 element
      int izero;
1059
      if ( basis_->real() )
Francois Gygi committed
1060 1061
        izero = 0;
      else
1062 1063 1064
        izero = basis_->rod_size(0)/2;
      //cout << " izero = " << izero << " G = " << basis_->kv(3*izero) << " " 
      //     << basis_->kv(3*izero+1) << " " << basis_->kv(3*izero+2) << endl;
Francois Gygi committed
1065 1066 1067
      p[izero] = complex<double> ( p[izero].real(), 0.0);
    }
    // reset values of empty rows of c_ to zero
1068
    for ( int i = basis_->localsize(); i < c_.mloc(); i++ )
Francois Gygi committed
1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086
    {
      p[i] = 0.0;
    }
  }
}

////////////////////////////////////////////////////////////////////////////////
SlaterDet& SlaterDet::operator=(SlaterDet& rhs)
{
  if ( this == &rhs ) return *this;
  assert(ctxt_.ictxt() == rhs.context().ictxt());
  c_ = rhs.c_;
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::memsize(void) const
{
1087
  return basis_->memsize() + c_.memsize();
Francois Gygi committed
1088 1089 1090 1091 1092
}

////////////////////////////////////////////////////////////////////////////////
double SlaterDet::localmemsize(void) const
{
1093
  return basis_->localmemsize() + c_.localmemsize();
Francois Gygi committed
1094 1095 1096 1097 1098
}

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::print(ostream& os, string encoding)
{
1099
  FourierTransform ft(*basis_,basis_->np(0),basis_->np(1),basis_->np(2));
Francois Gygi committed
1100 1101
  vector<complex<double> > wftmp(ft.np012loc());
  vector<double> wftmpr(ft.np012());
Francois Gygi committed
1102
  Base64Transcoder xcdr;
Francois Gygi committed
1103 1104 1105 1106 1107
  
  if ( ctxt_.onpe0() )
  {
    const double weight = 1.0; //!! fixed determinant weight to 1.0
    //!! no spin attribute written
1108
    os << "<slater_determinant kpoint=\"" << basis_->kpoint() << "\"\n"
Francois Gygi committed
1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128
       << "  weight=\"" << weight << "\""
       << " size=\"" << nst() << "\">" << endl;
 
    os << "<density_matrix form=\"diagonal\" size=\"" << nst() << "\">" 
       << 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;
  }
  
  for ( int n = 0; n < nst(); n++ )
  {
1129 1130 1131 1132
    // Barrier to limit the number of messages sent to task 0 
    // that don't have a receive posted
    ctxt_.barrier();
    
1133
    // transform data on ctxt_.mycol()
Francois Gygi committed
1134 1135 1136 1137 1138 1139 1140 1141 1142 1143
    if ( c_.pc(n) == ctxt_.mycol() )
    {
      //cout << " state " << n << " is stored on column " 
      //     << ctxt_.mycol() << " local index: " << c_.y(n) << endl;
      int nloc = c_.y(n); // local index
      ft.backward(c_.cvalptr(c_.mloc()*nloc),&wftmp[0]);
      
      double *a = (double*) &wftmp[0];
      for ( int i = 0; i < ft.np012loc(); i++ )
        wftmpr[i] = a[2*i];
1144 1145 1146 1147 1148 1149 1150 1151 1152 1153
    }
    
    // send blocks of wftmpr to pe0
    for ( int i = 0; i < ctxt_.nprow(); i++ )
    {
      bool iamsending = c_.pc(n) == ctxt_.mycol() && i == ctxt_.myrow();
      
      // send size of wftmpr block
      int size=-1;
      if ( ctxt_.onpe0() )
Francois Gygi committed
1154
      {
1155
        if ( iamsending )
Francois Gygi committed
1156
        {
1157 1158 1159 1160 1161 1162 1163 1164 1165 1166
          // sending to self, size not needed
        }
        else
          ctxt_.irecv(1,1,&size,1,i,c_.pc(n));
      }
      else
      {
        if ( iamsending )
        {
          size = ft.np012loc();
Francois Gygi committed
1167
          ctxt_.isend(1,1,&size,1,0,0);
1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187
        }
      }
      
      // 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);
          ctxt_.drecv(size,1,&wftmpr[istart],1,i,c_.pc(n));
        }
      }
      else
      {
        if ( iamsending )
        {
Francois Gygi committed
1188 1189 1190 1191
          ctxt_.dsend(size,1,&wftmpr[0],1,0,0);
        }
      }
    }
1192 1193

    // process the data      
Francois Gygi committed
1194 1195 1196 1197 1198 1199
    if ( ctxt_.onpe0() )
    {
      // wftmpr is now complete on task 0

      if ( encoding == "base64" )
      {
Francois Gygi committed
1200
        #if PLT_BIG_ENDIAN
Francois Gygi committed
1201
        xcdr.byteswap_double(ft.np012(),&wftmpr[0]);
Francois Gygi committed
1202
        #endif
1203 1204 1205
        int nbytes = ft.np012()*sizeof(double);
        int outlen = xcdr.nchars(nbytes);
        char* b = new char[outlen];
Francois Gygi committed
1206
        assert(b!=0);
1207
        xcdr.encode(nbytes,(byte*) &wftmpr[0],b);
Francois Gygi committed
1208 1209 1210 1211 1212
        // Note: optional x0,y0,z0 attributes not used, default is zero
        os << "<grid_function type=\"double\""
           << " nx=\"" << ft.np0()
           << "\" ny=\"" << ft.np1() << "\" nz=\"" << ft.np2() << "\""
           << " encoding=\"base64\">" << endl;
1213
        xcdr.print(outlen,(char*) b, os);
Francois Gygi committed
1214
        os << "</grid_function>\n";
Francois Gygi committed
1215
        delete [] b;
Francois Gygi committed
1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247
      }
      else
      {
        // encoding == "text" or unknown encoding
        // Note: optional x0,y0,z0 attributes not used, default is zero
        os << "<grid_function type=\"double\""
           << " nx=\"" << ft.np0()
           << "\" ny=\"" << ft.np1() << "\" nz=\"" << ft.np2() << "\""
           << " encoding=\"text\">" << endl;
        int count = 0;
        for ( int k = 0; k < ft.np2(); k++ )
          for ( int j = 0; j < ft.np1(); j++ )
            for ( int i = 0; i < ft.np0(); i++ )
            {
              os << " " << wftmpr[ft.index(i,j,k)];
              if ( count++%4 == 3)
                os << "\n";
            }
        if ( count%4 != 0 )
          os << "\n";
        os << "</grid_function>\n";
      }
    }
  }
  if ( ctxt_.onpe0() )
    os << "</slater_determinant>" << endl;
}

#if USE_CSTDIO_LFS
////////////////////////////////////////////////////////////////////////////////
void SlaterDet::write(FILE* outfile, string encoding)
{
1248
  FourierTransform ft(*basis_,basis_->np(0),basis_->np(1),basis_->np(2));
Francois Gygi committed
1249 1250
  vector<complex<double> > wftmp(ft.np012loc());
  vector<double> wftmpr(ft.np012());
Francois Gygi committed
1251
  Base64Transcoder xcdr;
Francois Gygi committed
1252 1253 1254 1255 1256 1257
  ostringstream os;
  
  if ( ctxt_.onpe0() )
  {
    const double weight = 1.0; //!! fixed determinant weight to 1.0
    //!! no spin attribute written
1258
    os << "<slater_determinant kpoint=\"" << basis_->kpoint() << "\"\n"
Francois Gygi committed
1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282
       << "  weight=\"" << weight << "\""
       << " size=\"" << nst() << "\">" << endl;
 
    os << "<density_matrix form=\"diagonal\" size=\"" << nst() << "\">" 
       << 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;
    
    string str(os.str());
    off_t len = str.length();
    fwrite(str.c_str(),sizeof(char),len,outfile);
  }
  
  for ( int n = 0; n < nst(); n++ )
  {
1283 1284 1285 1286
    // Barrier to limit the number of messages sent to task 0 
    // that don't have a receive posted
    ctxt_.barrier();
    
1287
    // transform data on ctxt_.mycol()
Francois Gygi committed
1288 1289 1290 1291 1292 1293 1294 1295 1296 1297
    if ( c_.pc(n) == ctxt_.mycol() )
    {
      //cout << " state " << n << " is stored on column " 
      //     << ctxt_.mycol() << " local index: " << c_.y(n) << endl;
      int nloc = c_.y(n); // local index
      ft.backward(c_.cvalptr(c_.mloc()*nloc),&wftmp[0]);
      
      double *a = (double*) &wftmp[0];
      for ( int i = 0; i < ft.np012loc(); i++ )
        wftmpr[i] = a[2*i];
1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316
    }
    
    // send blocks of wftmpr to pe0
    for ( int i = 0; i < ctxt_.nprow(); i++ )
    {
      bool iamsending = c_.pc(n) == ctxt_.mycol() && i == ctxt_.myrow();
      
      // send size of wftmpr block
      int size=-1;
      if ( ctxt_.onpe0() )
      {
        if ( iamsending )
        {
          // sending to self, size not needed
        }
        else
          ctxt_.irecv(1,1,&size,1,i,c_.pc(n));
      }
      else
Francois Gygi committed
1317
      {
1318
        if ( iamsending )
Francois Gygi committed
1319
        {
1320
          size = ft.np012loc();
Francois Gygi committed
1321
          ctxt_.isend(1,1,&size,1,0,0);
1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341
        }
      }
      
      // 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);
          ctxt_.drecv(size,1,&wftmpr[istart],1,i,c_.pc(n));
        }
      }
      else
      {
        if ( iamsending )
        {
Francois Gygi committed
1342 1343 1344 1345
          ctxt_.dsend(size,1,&wftmpr[0],1,0,0);
        }
      }
    }
1346 1347
      
    // process data
Francois Gygi committed
1348 1349 1350 1351 1352 1353
    if ( ctxt_.onpe0() )
    {
      // wftmpr is now complete on task 0

      if ( encoding == "base64" )
      {
Francois Gygi committed
1354
        #if PLT_BIG_ENDIAN
Francois Gygi committed
1355
        xcdr.byteswap_double(ft.np012(),&wftmpr[0]);
Francois Gygi committed
1356
        #endif
1357 1358 1359 1360
        
        int nbytes = ft.np012()*sizeof(double);
        int outlen = xcdr.nchars(nbytes);
        char* b = new char[outlen];
Francois Gygi committed
1361
        assert(b!=0);
1362
        xcdr.encode(nbytes,(byte*) &wftmpr[0],b);
Francois Gygi committed
1363
        
Francois Gygi committed
1364 1365 1366 1367 1368 1369
        // Note: optional x0,y0,z0 attributes not used, default is zero
        os.str("");
        os << "<grid_function type=\"double\""
           << " nx=\"" << ft.np0()
           << "\" ny=\"" << ft.np1() << "\" nz=\"" << ft.np2() << "\""
           << " encoding=\"base64\">" << endl;
1370
           
Francois Gygi committed
1371 1372 1373
        string str(os.str());
        off_t len = str.length();
        fwrite(str.c_str(),sizeof(char),len,outfile);
1374 1375 1376 1377
        
        // write buffer b inserting newlines 
        xcdr.print(outlen, b, outfile);
        
Francois Gygi committed
1378 1379 1380 1381 1382
        os.str("");
        os << "</grid_function>\n";
        str = os.str();
        len = str.length();
        fwrite(str.c_str(),sizeof(char),len,outfile);
Francois Gygi committed
1383
        delete [] b;
Francois Gygi committed
1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443
      }
      else
      {
        // encoding == "text" or unknown encoding
        // Note: optional x0,y0,z0 attributes not used, default is zero
        os.str("");
        os << "<grid_function type=\"double\""
           << " nx=\"" << ft.np0()
           << "\" ny=\"" << ft.np1() << "\" nz=\"" << ft.np2() << "\""
           << " encoding=\"text\">" << endl;
        string str(os.str());
        off_t len = str.length();
        fwrite(str.c_str(),sizeof(char),len,outfile);
        int count = 0;
        for ( int k = 0; k < ft.np2(); k++ )
        {
          os.str("");
          for ( int j = 0; j < ft.np1(); j++ )
          {
            for ( int i = 0; i < ft.np0(); i++ )
            {
              os << " " << wftmpr[ft.index(i,j,k)];
              if ( count++%4 == 3)
                os << "\n";
            }
          }
          str = os.str();
          len = str.length();
          fwrite(str.c_str(),sizeof(char),len,outfile);
        }
        os.str("");
        if ( count%4 != 0 )
          os << "\n";
        os << "</grid_function>\n";
        str = os.str();
        len = str.length();
        fwrite(str.c_str(),sizeof(char),len,outfile);
      }

    }
  }
  
  if ( ctxt_.onpe0() )
  {
    os.str("");
    os << "</slater_determinant>" << endl;
    string str(os.str());
    off_t len = str.length();
    fwrite(str.c_str(),sizeof(char),len,outfile);
  }
}
#endif

////////////////////////////////////////////////////////////////////////////////
void SlaterDet::info(ostream& os)
{  
  if ( ctxt_.onpe0() )
  {
    const double weight = 1.0; //!! fixed determinant weight to 1.0
    //!! no spin attribute written
1444
    os << "<slater_determinant kpoint=\"" << basis_->kpoint() << "\"\n"
Francois Gygi committed
1445 1446
       << "  weight=\"" << weight << "\""
       << " size=\"" << nst() << "\">" << endl;
1447 1448
    os << " <!-- sdcontext: " << ctxt_.nprow() << "x" << ctxt_.npcol() << " -->"
       << endl;
1449 1450
    //os << " <!-- sdcontext: " << ctxt_ << " -->" << endl;
    os << " <!-- basis size: " << basis_->size() << " -->" << endl;
1451
    os << " <!-- c dimensions: "
Francois Gygi committed
1452 1453
       << c_.m() << "x" << c_.n()
       << "   (" << c_.mb() << "x" << c_.nb() << " blocks)" << " -->" << endl;
1454
    os << " <density_matrix form=\"diagonal\" size=\"" << nst() << "\">" 
Francois Gygi committed
1455
       << endl;
1456
    os << " </density_matrix>" << endl;
Francois Gygi committed
1457 1458 1459 1460 1461 1462 1463 1464 1465 1466
    os << "</slater_determinant>" << endl;
  }
}

////////////////////////////////////////////////////////////////////////////////
ostream& operator<<(ostream& os, SlaterDet& sd)
{
  sd.print(os,"text");
  return os;
}