NonLocalPotential.C 60.8 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
// NonLocalPotential.C
//
////////////////////////////////////////////////////////////////////////////////
18
// $Id: NonLocalPotential.C,v 1.28 2009-11-30 02:32:13 fgygi Exp $
Francois Gygi committed
19 20

#include "NonLocalPotential.h"
Francois Gygi committed
21
#include "Species.h"
Francois Gygi committed
22
#include "blas.h"
23 24
#include <iomanip>
using namespace std;
Francois Gygi committed
25

26
#if USE_MASSV
Francois Gygi committed
27 28 29 30 31 32 33
extern "C" void vsincos(double *x, double *y, double *z, int *n);
#endif


////////////////////////////////////////////////////////////////////////////////
NonLocalPotential::~NonLocalPotential(void)
{
Francois Gygi committed
34
#if TIMING
35 36 37 38 39 40 41 42 43
  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 )
    {
44 45 46 47 48
      cout << "<timing name=\""
           << setw(15) << (*i).first << "\""
           << " min=\"" << setprecision(3) << setw(9) << tmin << "\""
           << " max=\"" << setprecision(3) << setw(9) << tmax << "\"/>"
           << endl;
49 50
    }
  }
Francois Gygi committed
51
#endif
Francois Gygi committed
52 53 54 55 56 57
}

////////////////////////////////////////////////////////////////////////////////
void NonLocalPotential::init(void)
{
  const int ngwl = basis_.localsize();
58

Francois Gygi committed
59
  nsp = atoms_.nsp();
60

Francois Gygi committed
61 62 63 64 65 66 67 68 69
  lmax.resize(nsp);
  lloc.resize(nsp);
  lproj.resize(nsp);
  na.resize(nsp);
  npr.resize(nsp);
  nprna.resize(nsp);
  wt.resize(nsp);
  twnl.resize(nsp);
  dtwnl.resize(nsp);
70

Francois Gygi committed
71 72 73
  nquad.resize(nsp);
  rquad.resize(nsp);
  wquad.resize(nsp);
74

Francois Gygi committed
75 76 77 78
  nspnl = 0;
  for ( int is = 0; is < nsp; is++ )
  {
    Species *s = atoms_.species_list[is];
79

Francois Gygi committed
80 81
    npr[is] = 0;
    nprna[is] = 0;
82

Francois Gygi committed
83 84 85 86 87 88 89
    if ( s->non_local() )
    {
      nspnl++;
      na[is] = atoms_.na(is);
      lmax[is] = s->lmax();
      lloc[is] = s->llocal();
      nquad[is] = s->nquad();
90

Francois Gygi committed
91 92 93 94 95 96 97 98 99 100 101 102
      // compute total number of projectors npr[is]
      // KB potentials have nlm projectors
      // semilocal potentials have nlm*nquad projectors
      if ( s->nquad() == 0 )
      {
        npr[is] = s->nlm();
      }
      else
      {
        npr[is] = s->nlm() * nquad[is];
      }
      nprna[is] = npr[is] * na[is];
103

Francois Gygi committed
104 105
      // l value for projector ipr
      lproj[is].resize(npr[is]);
106 107

      twnl[is].resize(npr[is]*ngwl);
108
      dtwnl[is].resize(npr[is]*6*ngwl);
109

Francois Gygi committed
110
      // quadrature abcissae and weights
Francois Gygi committed
111 112
      rquad[is].resize(nquad[is]);
      wquad[is].resize(nquad[is]);
113 114 115

      enum quadrature_rule_type { TRAPEZOID, MODIF_TRAPEZOID,
      TRAPEZOID_WITH_RIGHT_ENDPOINT,
116
      SIMPSON };
117

118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161
      const quadrature_rule_type quad_rule = TRAPEZOID;
      //const quadrature_rule_type quad_rule = MODIF_TRAPEZOID;
      //const quadrature_rule_type quad_rule = TRAPEZOID_WITH_RIGHT_ENDPOINT;
      //const quadrature_rule_type quad_rule = SIMPSON;

      if ( quad_rule == TRAPEZOID )
      {
        // trapezoidal rule with interior points only
        // (end points are zero)
        const double h = s->rquad() / (nquad[is]+1);
        for ( int iquad = 0; iquad < nquad[is]; iquad++ )
        {
          rquad[is][iquad] = (iquad+1) * h;
          wquad[is][iquad] = h;
        }
        //cout << " NonLocalPotential::init: trapezoidal rule (interior)"
        //     << endl;
      }
      else if ( quad_rule == MODIF_TRAPEZOID )
      {
        // use modified trapezoidal rule with interior points, and include
        // correction for first derivative at r=0 as
        // h^2/12 f'(0) where f'(0) is estimated with f(h)/h
        // i.e. add correction h/12) * f(h)
        // See Davis & Rabinowitz, p. 132
        const double h = s->rquad() / (nquad[is]+1);
        for ( int iquad = 0; iquad < nquad[is]; iquad++ )
        {
          rquad[is][iquad] = (iquad+1) * h;
          wquad[is][iquad] = h;
        }
        wquad[is][0] += h / 12.0;
        //cout << " NonLocalPotential::init: modified trapezoidal rule"
        //     << endl;
      }
      else if ( quad_rule == TRAPEZOID_WITH_RIGHT_ENDPOINT )
      {
        const double h = s->rquad() / nquad[is];
        for ( int iquad = 0; iquad < nquad[is]; iquad++ )
        {
          rquad[is][iquad] = (iquad+1) * h;
          wquad[is][iquad] = h;
        }
        wquad[is][nquad[is]-1] = 0.5 * h;
Francois Gygi committed
162
        //cout << " NonLocalPotential::init: trapezoidal rule, right endpoint"
163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182
        //     << endl;
      }
      else if ( quad_rule == SIMPSON )
      {
        // must have 2n+1 points
        assert(nquad[is]%2==1);
        const double h = s->rquad() / (nquad[is]-1);
        for ( int iquad = 0; iquad < nquad[is]; iquad++ )
        {
          rquad[is][iquad] = iquad * h;
          if ( ( iquad == 0 ) || ( iquad == nquad[is]-1 ) )
            wquad[is][iquad] = h / 3.0;
          else if ( iquad % 2 == 0 )
            wquad[is][iquad] = h * 2.0 / 3.0;
          else
            wquad[is][iquad] = h * 4.0 / 3.0;
        }
        //cout << " NonLocalPotential::init: Simpson rule" << endl;
      }
      else
Francois Gygi committed
183
      {
184
        assert(false);
Francois Gygi committed
185
      }
186

Francois Gygi committed
187 188
      // compute weights wt[is][ipr]
      wt[is].resize(npr[is]);
189 190

      // compute lproj[is][ipr]
Francois Gygi committed
191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214
      int ipr_base = 0;
      for ( int l = 0; l <= lmax[is]; l++ )
      {
        if ( l != lloc[is] )
        {
          if ( nquad[is] == 0 )
          {
            // Kleinman-Bylander form
            // wt[is][ipr]
            // index = ipr_base+m
            for ( int m = 0; m < 2*l+1; m++ )
            {
              const int ipr = ipr_base + m;
              wt[is][ipr] = s->wsg(l);
              lproj[is][ipr] = l;
            }
            ipr_base += 2*l+1;
          }
          else
          {
            for ( int iquad = 0; iquad < nquad[is]; iquad++ )
            {
              const double r = rquad[is][iquad];
              double v,dv,vl,dvl;
215 216
              s->dvpsr(l,r,v,dv);
              s->dvpsr(lloc[is],r,vl,dvl);
Francois Gygi committed
217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240
              // wt[is][iquad+ipr*nquad]
              for ( int m = 0; m < 2*l+1; m++ )
              {
                const int ipr = ipr_base + iquad + nquad[is] * m;
                wt[is][ipr] = ( v - vl ) * wquad[is][iquad];
                lproj[is][ipr] = l;
              }
            }
            ipr_base += (2*l+1) * nquad[is];
          }
        }
      }
      assert(ipr_base==npr[is]);
    } // if s->non_local()
  }
}

////////////////////////////////////////////////////////////////////////////////
void NonLocalPotential::update_twnl(void)
{
  // update arrays twnl[is][ipr][ig], dtwnl[is][ipr][j][ig],
  // following a change of cell dimensions
  // It is assumed that basis_ has been updated
  // It is assumed that nsp, npr[is], nquad[is] did not change since init
241

242
  tmap["update_twnl"].start();
243

Francois Gygi committed
244 245 246 247
  const int ngwl = basis_.localsize();
  const double pi = M_PI;
  const double fpi = 4.0 * pi;
  const double s14pi = sqrt(1.0/fpi);
248
  const double s34pi = sqrt(3.0/fpi);
Francois Gygi committed
249
  const double s54pi = sqrt(5.0/fpi);
250 251
  const double s20pi = sqrt(20.0*pi);
  const double s20pi3 = sqrt(20.0*pi/3.0);
Francois Gygi committed
252 253
  const double s3 = sqrt(3.0);

Francois Gygi committed
254 255 256 257 258 259
  const double *kpg   = basis_.kpg_ptr();
  const double *kpg2  = basis_.kpg2_ptr();
  const double *kpgi  = basis_.kpgi_ptr();
  const double *kpg_x = basis_.kpgx_ptr(0);
  const double *kpg_y = basis_.kpgx_ptr(1);
  const double *kpg_z = basis_.kpgx_ptr(2);
260

Francois Gygi committed
261 262 263 264
  // compute twnl and dtwnl
  for ( int is = 0; is < nsp; is++ )
  {
    Species *s = atoms_.species_list[is];
265

Francois Gygi committed
266 267
    int ilm = 0;
    for ( int l = 0; l <= lmax[is]; l++ )
268
    {
Francois Gygi committed
269 270 271 272 273 274 275
      if ( l != lloc[is] )
      {
        if ( l == 0 )
        {
          if ( nquad[is] == 0 )
          {
            // Kleinman-Bylander
276

Francois Gygi committed
277 278 279 280
            // twnl[is][ipr][ig]
            // ipr = ilm = 0
            // index = ig + ngwl*ipr, i.e. index = ig
            double *t0  = &twnl[is][0];
281

282 283 284 285 286 287 288 289 290
            // dtwnl[is][ipr][ij][ngwl]
            // index = ig + ngwl * ( ij + 6 * ipr ), ipr = 0
            // i.e. index = ig + ij * ngwl
            double *dt0_xx = &dtwnl[is][0*ngwl];
            double *dt0_yy = &dtwnl[is][1*ngwl];
            double *dt0_zz = &dtwnl[is][2*ngwl];
            double *dt0_xy = &dtwnl[is][3*ngwl];
            double *dt0_yz = &dtwnl[is][4*ngwl];
            double *dt0_xz = &dtwnl[is][5*ngwl];
Francois Gygi committed
291
            // Special case k=G=0 is ok since kpgi[0] = 0.0 at k=G=0
Francois Gygi committed
292 293 294
            for ( int ig = 0; ig < ngwl; ig++ )
            {
              double v,dv;
Francois Gygi committed
295
              s->dvnlg(0,kpg[ig],v,dv);
Francois Gygi committed
296 297

              t0[ig] = s14pi * v;
298

Francois Gygi committed
299 300 301
              const double tgx = kpg_x[ig];
              const double tgy = kpg_y[ig];
              const double tgz = kpg_z[ig];
Francois Gygi committed
302 303 304
              const double tgx2 = tgx * tgx;
              const double tgy2 = tgy * tgy;
              const double tgz2 = tgz * tgz;
305

Francois Gygi committed
306
              const double tmp = kpgi[ig] * s14pi * dv;
307 308 309 310 311 312
              dt0_xx[ig] = tmp * tgx * tgx;
              dt0_yy[ig] = tmp * tgy * tgy;
              dt0_zz[ig] = tmp * tgz * tgz;
              dt0_xy[ig] = tmp * tgx * tgy;
              dt0_yz[ig] = tmp * tgy * tgz;
              dt0_xz[ig] = tmp * tgx * tgz;
Francois Gygi committed
313 314 315 316 317 318 319 320 321 322 323 324 325
            }
          }
          else
          {
            // semi-local
            for ( int iquad = 0; iquad < nquad[is]; iquad++ )
            {
              // twnl[is][ipr][ig]
              // ipr = iquad + nquad[is]*ilm, where ilm=0
              //     = iquad
              // index = ig + ngwl*iquad
              double *t0 = &twnl[is][ngwl*iquad];
              // dtwnl[is][ipr][j][ngwl]
326
              // index = ig + ngwl * ( ij + 6 * iquad)
327 328 329 330 331 332
              double *dt0_xx = &dtwnl[is][ngwl*(0+6*iquad)];
              double *dt0_yy = &dtwnl[is][ngwl*(1+6*iquad)];
              double *dt0_zz = &dtwnl[is][ngwl*(2+6*iquad)];
              double *dt0_xy = &dtwnl[is][ngwl*(3+6*iquad)];
              double *dt0_yz = &dtwnl[is][ngwl*(4+6*iquad)];
              double *dt0_xz = &dtwnl[is][ngwl*(5+6*iquad)];
Francois Gygi committed
333
              const double r = rquad[is][iquad];
Francois Gygi committed
334 335 336

              // G=0 element must be treated separately if k=0
              if ( basis_.real() && ctxt_.myrow() == 0)
337
              {
Francois Gygi committed
338 339
                // compute G=0 element separately to get
                // sin(Gr)/(Gr) limit -> 1
340
                const int ig = 0;
341
                // this task holds the G=0 element
342 343
                // I(l=0) = 4 pi j_l(G r) r
                // twnl[is][ipr][l][ig] = 4 pi j_0(Gr_i) r_i Ylm
344
                // j_0(Gr) = sin(Gr) / (Gr) == 1.0
345
                // Ylm = s14pi
346

347 348
                // next line: sin(Gr)/(Gr) replaced by 1
                t0[ig] = fpi * s14pi * r;
349

350
                // dtwnl = fpi s14pi G_i G_j / G (r cos(Gr)/G -sin(Gr)/G^2)
351 352 353 354 355 356 357
                // dtwnl = 0.0;
                dt0_xx[ig] = 0.0;
                dt0_yy[ig] = 0.0;
                dt0_zz[ig] = 0.0;
                dt0_xy[ig] = 0.0;
                dt0_yz[ig] = 0.0;
                dt0_xz[ig] = 0.0;
358 359 360
              }
              else
              {
Francois Gygi committed
361
                // Use the normal procedure
362 363 364 365 366
                const int ig = 0;
                // I(l=0) = 4 pi j_l(G r) r
                // twnl[is][ipr][l][ig] = 4 pi j_0(Gr_i) r_i Ylm
                // j_0(Gr) * r = sin(Gr) / G
                // Ylm = s14pi
Francois Gygi committed
367
                const double arg = kpg[ig] * r;
368
                // Note: for G=0, gi[0] = 0
369

Francois Gygi committed
370 371 372 373
                const double tgx = kpg_x[ig];
                const double tgy = kpg_y[ig];
                const double tgz = kpg_z[ig];
                const double tgi = kpgi[ig];
374
                const double tgi2 = tgi * tgi;
375

376 377
                const double ts = sin(arg);
                const double tc = cos(arg);
378

379
                t0[ig] = fpi * s14pi * ts * tgi;
380

Francois Gygi committed
381 382
                // dtwnl = fpi s14pi k+G_i k+G_j / k+G
                // (r cos(k+Gr)/k+G -sin(k+Gr)/k+G^2)
383 384 385 386 387 388 389 390 391 392
                const double tmp = fpi * s14pi * tgi2 * (r*tc - ts*tgi);
                dt0_xx[ig] = tmp * tgx * tgx;
                dt0_yy[ig] = tmp * tgy * tgy;
                dt0_zz[ig] = tmp * tgz * tgz;
                dt0_xy[ig] = tmp * tgx * tgy;
                dt0_yz[ig] = tmp * tgy * tgz;
                dt0_xz[ig] = tmp * tgx * tgz;
              }

              for ( int ig = 1; ig < ngwl; ig++ )
Francois Gygi committed
393
              {
Francois Gygi committed
394 395 396
                // I(l=0) = 4 pi j_l(|k+G| r) r
                // twnl[is][ipr][l][ig] = 4 pi j_0(|k+G|r_i) r_i Ylm
                // j_0(|k+G|r) * r = sin(|k+G|r) / |k+G|
Francois Gygi committed
397
                // Ylm = s14pi
Francois Gygi committed
398
                const double arg = kpg[ig] * r;
399

Francois Gygi committed
400 401 402 403
                const double tgx = kpg_x[ig];
                const double tgy = kpg_y[ig];
                const double tgz = kpg_z[ig];
                const double tgi = kpgi[ig];
Francois Gygi committed
404
                const double tgi2 = tgi * tgi;
405

Francois Gygi committed
406 407
                const double ts = sin(arg);
                const double tc = cos(arg);
408

Francois Gygi committed
409
                t0[ig] = fpi * s14pi * ts * tgi;
410

Francois Gygi committed
411 412
                // dtwnl = fpi s14pi (k+G)_i (k+G)_j /
                // |k+G| (r cos(|k+G|r)/|k+G| -sin(|k+G|r)/|k+G|^2)
413 414 415 416 417 418 419
                const double tmp = fpi * s14pi * tgi2 * (r*tc - ts*tgi);
                dt0_xx[ig] = tmp * tgx * tgx;
                dt0_yy[ig] = tmp * tgy * tgy;
                dt0_zz[ig] = tmp * tgz * tgz;
                dt0_xy[ig] = tmp * tgx * tgy;
                dt0_yz[ig] = tmp * tgy * tgz;
                dt0_xz[ig] = tmp * tgx * tgz;
Francois Gygi committed
420
              }
Francois Gygi committed
421
            } // for iquad
Francois Gygi committed
422 423 424 425 426 427 428 429
          }
          ilm += 2*l+1;
        }
        else if ( l == 1 )
        {
          if ( nquad[is] == 0 )
          {
            // Kleinman-Bylander
430

Francois Gygi committed
431 432 433 434 435 436 437 438 439
            // twnl[is][ipr][ig]
            // ipr = ilm
            const int ipr1 = ilm;
            const int ipr2 = ilm+1;
            const int ipr3 = ilm+2;
            // index = ig + ngwl*ilm
            double *t1 = &twnl[is][ngwl*ipr1];
            double *t2 = &twnl[is][ngwl*ipr2];
            double *t3 = &twnl[is][ngwl*ipr3];
440

441 442 443 444 445 446 447 448 449 450 451 452 453 454 455
            // dtwnl[is][ipr][ij][ngwl]
            // index = ig + ngwl * ( ij + 6 * ipr )
            double *dt1_xx = &dtwnl[is][ngwl*(0+6*ipr1)];
            double *dt1_yy = &dtwnl[is][ngwl*(1+6*ipr1)];
            double *dt1_zz = &dtwnl[is][ngwl*(2+6*ipr1)];
            double *dt1_xy = &dtwnl[is][ngwl*(3+6*ipr1)];
            double *dt1_yz = &dtwnl[is][ngwl*(4+6*ipr1)];
            double *dt1_xz = &dtwnl[is][ngwl*(5+6*ipr1)];

            double *dt2_xx = &dtwnl[is][ngwl*(0+6*ipr2)];
            double *dt2_yy = &dtwnl[is][ngwl*(1+6*ipr2)];
            double *dt2_zz = &dtwnl[is][ngwl*(2+6*ipr2)];
            double *dt2_xy = &dtwnl[is][ngwl*(3+6*ipr2)];
            double *dt2_yz = &dtwnl[is][ngwl*(4+6*ipr2)];
            double *dt2_xz = &dtwnl[is][ngwl*(5+6*ipr2)];
Francois Gygi committed
456

457 458 459 460 461 462
            double *dt3_xx = &dtwnl[is][ngwl*(0+6*ipr3)];
            double *dt3_yy = &dtwnl[is][ngwl*(1+6*ipr3)];
            double *dt3_zz = &dtwnl[is][ngwl*(2+6*ipr3)];
            double *dt3_xy = &dtwnl[is][ngwl*(3+6*ipr3)];
            double *dt3_yz = &dtwnl[is][ngwl*(4+6*ipr3)];
            double *dt3_xz = &dtwnl[is][ngwl*(5+6*ipr3)];
Francois Gygi committed
463 464 465 466

            for ( int ig = 0; ig < ngwl; ig++ )
            {
              double v,dv;
Francois Gygi committed
467
              const double tg = kpg[ig];
468
              s->dvnlg(l,tg,v,dv);
469

Francois Gygi committed
470 471 472
              const double tgx = kpg_x[ig];
              const double tgy = kpg_y[ig];
              const double tgz = kpg_z[ig];
Francois Gygi committed
473 474 475
              const double tgx2 = tgx * tgx;
              const double tgy2 = tgy * tgy;
              const double tgz2 = tgz * tgz;
476

Francois Gygi committed
477
              const double tgi = kpgi[ig];
Francois Gygi committed
478
              const double tgi2 = tgi * tgi;
479

Francois Gygi committed
480 481 482
              const double y1 = s34pi * tgx * tgi;
              const double y2 = s34pi * tgy * tgi;
              const double y3 = s34pi * tgz * tgi;
483

Francois Gygi committed
484 485 486 487
              t1[ig]  = y1 * v;
              t2[ig]  = y2 * v;
              t3[ig]  = y3 * v;

488 489 490 491 492 493 494 495
              const double fac1 = - y1 * ( v - tg * dv ) * tgi2;
              // m=x
              dt1_xx[ig] = fac1 * tgx2 + v * y1;
              dt1_yy[ig] = fac1 * tgy2;
              dt1_zz[ig] = fac1 * tgz2;
              dt1_xy[ig] = fac1 * tgx * tgy;
              dt1_yz[ig] = fac1 * tgy * tgz;
              dt1_xz[ig] = fac1 * tgx * tgz;
496

497 498 499 500 501 502 503 504
              const double fac2 = - y2 * ( v - tg * dv ) * tgi2;
              // m=y
              dt2_xx[ig] = fac2 * tgx2;
              dt2_yy[ig] = fac2 * tgy2 + v * y2;
              dt2_zz[ig] = fac2 * tgz2;
              dt2_xy[ig] = fac2 * tgx * tgy + v * y1;
              dt2_yz[ig] = fac2 * tgy * tgz;
              dt2_xz[ig] = fac2 * tgx * tgz;
505

506 507 508 509 510 511 512 513
              const double fac3 = - y3 * ( v - tg * dv ) * tgi2;
              // m=z
              dt3_xx[ig] = fac3 * tgx2;
              dt3_yy[ig] = fac3 * tgy2;
              dt3_zz[ig] = fac3 * tgz2 + v * y3;
              dt3_xy[ig] = fac3 * tgx * tgy;
              dt3_yz[ig] = fac3 * tgy * tgz + v * y2;
              dt3_xz[ig] = fac3 * tgx * tgz + v * y1;
Francois Gygi committed
514 515 516 517 518 519 520 521 522 523 524 525 526 527 528
            }
          }
          else
          {
            // semi-local
            for ( int iquad = 0; iquad < nquad[is]; iquad++ )
            {
              // twnl[is][ipr][ig]
              // index = ig + ngwl*(iquad+nquad[is]*ilm)
              const int ipr1 = iquad+nquad[is]*ilm;
              const int ipr2 = iquad+nquad[is]*(ilm+1);
              const int ipr3 = iquad+nquad[is]*(ilm+2);
              double *t1 = &twnl[is][ngwl*ipr1];
              double *t2 = &twnl[is][ngwl*ipr2];
              double *t3 = &twnl[is][ngwl*ipr3];
529

Francois Gygi committed
530
              // dtwnl[is][ipr][j][ngwl]
531 532 533 534 535 536 537
              // index = ig + ngwl * ( ij + 6 * ipr )
              double *dt1_xx = &dtwnl[is][ngwl*(0+6*ipr1)];
              double *dt1_yy = &dtwnl[is][ngwl*(1+6*ipr1)];
              double *dt1_zz = &dtwnl[is][ngwl*(2+6*ipr1)];
              double *dt1_xy = &dtwnl[is][ngwl*(3+6*ipr1)];
              double *dt1_yz = &dtwnl[is][ngwl*(4+6*ipr1)];
              double *dt1_xz = &dtwnl[is][ngwl*(5+6*ipr1)];
Francois Gygi committed
538

539 540 541 542 543 544
              double *dt2_xx = &dtwnl[is][ngwl*(0+6*ipr2)];
              double *dt2_yy = &dtwnl[is][ngwl*(1+6*ipr2)];
              double *dt2_zz = &dtwnl[is][ngwl*(2+6*ipr2)];
              double *dt2_xy = &dtwnl[is][ngwl*(3+6*ipr2)];
              double *dt2_yz = &dtwnl[is][ngwl*(4+6*ipr2)];
              double *dt2_xz = &dtwnl[is][ngwl*(5+6*ipr2)];
Francois Gygi committed
545

546 547 548 549 550 551
              double *dt3_xx = &dtwnl[is][ngwl*(0+6*ipr3)];
              double *dt3_yy = &dtwnl[is][ngwl*(1+6*ipr3)];
              double *dt3_zz = &dtwnl[is][ngwl*(2+6*ipr3)];
              double *dt3_xy = &dtwnl[is][ngwl*(3+6*ipr3)];
              double *dt3_yz = &dtwnl[is][ngwl*(4+6*ipr3)];
              double *dt3_xz = &dtwnl[is][ngwl*(5+6*ipr3)];
552

Francois Gygi committed
553 554 555 556 557
              const double r = rquad[is][iquad];
              for ( int ig = 0; ig < ngwl; ig++ )
              {
                double v = 0.0, dv = 0.0;
                // j_1(Gr) = (1/(Gr))*(sin(Gr)/(Gr)-cos(Gr))
Francois Gygi committed
558
                const double tg = kpg[ig];
559
                const double z = tg * r;
Francois Gygi committed
560 561 562 563 564 565
                if ( z != 0.0 )
                {
                  const double zi = 1.0 / z;
                  const double c = cos(z);
                  const double s = sin(z);
                  const double j1 = ( s * zi - c ) * zi;
566
                  const double dj1 =
Francois Gygi committed
567 568 569 570 571 572 573
                    ( 2.0 * z * c + ( z*z - 2.0 ) * s ) * zi*zi*zi;
                  // v = 4 pi j1(Gr) r
                  v = fpi * j1 * r;
                  // dv = d/dG v = 4 pi dj1(Gr)/d(Gr) d(Gr)/dG r
                  //    = 4 pi dj1 r^2
                  dv = fpi * dj1 * r * r;
                }
574

Francois Gygi committed
575 576 577
                const double tgx = kpg_x[ig];
                const double tgy = kpg_y[ig];
                const double tgz = kpg_z[ig];
Francois Gygi committed
578 579 580
                const double tgx2 = tgx * tgx;
                const double tgy2 = tgy * tgy;
                const double tgz2 = tgz * tgz;
581

Francois Gygi committed
582
                const double tgi = kpgi[ig];
Francois Gygi committed
583
                const double tgi2 = tgi * tgi;
584

Francois Gygi committed
585 586 587
                const double y1 = s34pi * tgx * tgi;
                const double y2 = s34pi * tgy * tgi;
                const double y3 = s34pi * tgz * tgi;
588

Francois Gygi committed
589 590 591 592
                t1[ig]  = y1 * v;
                t2[ig]  = y2 * v;
                t3[ig]  = y3 * v;

593 594 595 596 597 598 599 600
                const double fac1 = - y1 * ( v - tg * dv ) * tgi2;
                // m=x
                dt1_xx[ig] = fac1 * tgx2 + v * y1;
                dt1_yy[ig] = fac1 * tgy2;
                dt1_zz[ig] = fac1 * tgz2;
                dt1_xy[ig] = fac1 * tgx * tgy;
                dt1_yz[ig] = fac1 * tgy * tgz;
                dt1_xz[ig] = fac1 * tgx * tgz;
601

602 603 604 605 606 607 608 609
                const double fac2 = - y2 * ( v - tg * dv ) * tgi2;
                // m=y
                dt2_xx[ig] = fac2 * tgx2;
                dt2_yy[ig] = fac2 * tgy2 + v * y2;
                dt2_zz[ig] = fac2 * tgz2;
                dt2_xy[ig] = fac2 * tgx * tgy + v * y1;
                dt2_yz[ig] = fac2 * tgy * tgz;
                dt2_xz[ig] = fac2 * tgx * tgz;
610

611 612 613 614 615 616 617 618
                const double fac3 = - y3 * ( v - tg * dv ) * tgi2;
                // m=z
                dt3_xx[ig] = fac3 * tgx2;
                dt3_yy[ig] = fac3 * tgy2;
                dt3_zz[ig] = fac3 * tgz2 + v * y3;
                dt3_xy[ig] = fac3 * tgx * tgy;
                dt3_yz[ig] = fac3 * tgy * tgz + v * y2;
                dt3_xz[ig] = fac3 * tgx * tgz + v * y1;
Francois Gygi committed
619 620 621 622 623 624 625 626 627 628 629 630 631 632 633
              } // ig
            }
          }
          ilm += 2*l+1;
        }
        else if ( l == 2 )
        {
          if ( nquad[is] == 0 )
          {
            // Kleinman-Bylander
            const int ipr4 = ilm;
            const int ipr5 = ilm+1;
            const int ipr6 = ilm+2;
            const int ipr7 = ilm+3;
            const int ipr8 = ilm+4;
634

Francois Gygi committed
635 636 637 638 639
            double *t4 = &twnl[is][ngwl*ipr4];
            double *t5 = &twnl[is][ngwl*ipr5];
            double *t6 = &twnl[is][ngwl*ipr6];
            double *t7 = &twnl[is][ngwl*ipr7];
            double *t8 = &twnl[is][ngwl*ipr8];
640

641 642 643 644 645 646 647 648
            // dtwnl[is][ipr][ij][ngwl]
            // index = ig + ngwl * ( ij + 6 * ipr )
            double *dt4_xx = &dtwnl[is][ngwl*(0+6*ipr4)];
            double *dt4_yy = &dtwnl[is][ngwl*(1+6*ipr4)];
            double *dt4_zz = &dtwnl[is][ngwl*(2+6*ipr4)];
            double *dt4_xy = &dtwnl[is][ngwl*(3+6*ipr4)];
            double *dt4_yz = &dtwnl[is][ngwl*(4+6*ipr4)];
            double *dt4_xz = &dtwnl[is][ngwl*(5+6*ipr4)];
Francois Gygi committed
649

650 651 652 653 654 655
            double *dt5_xx = &dtwnl[is][ngwl*(0+6*ipr5)];
            double *dt5_yy = &dtwnl[is][ngwl*(1+6*ipr5)];
            double *dt5_zz = &dtwnl[is][ngwl*(2+6*ipr5)];
            double *dt5_xy = &dtwnl[is][ngwl*(3+6*ipr5)];
            double *dt5_yz = &dtwnl[is][ngwl*(4+6*ipr5)];
            double *dt5_xz = &dtwnl[is][ngwl*(5+6*ipr5)];
656

657 658 659 660 661 662
            double *dt6_xx = &dtwnl[is][ngwl*(0+6*ipr6)];
            double *dt6_yy = &dtwnl[is][ngwl*(1+6*ipr6)];
            double *dt6_zz = &dtwnl[is][ngwl*(2+6*ipr6)];
            double *dt6_xy = &dtwnl[is][ngwl*(3+6*ipr6)];
            double *dt6_yz = &dtwnl[is][ngwl*(4+6*ipr6)];
            double *dt6_xz = &dtwnl[is][ngwl*(5+6*ipr6)];
663

664 665 666 667 668 669
            double *dt7_xx = &dtwnl[is][ngwl*(0+6*ipr7)];
            double *dt7_yy = &dtwnl[is][ngwl*(1+6*ipr7)];
            double *dt7_zz = &dtwnl[is][ngwl*(2+6*ipr7)];
            double *dt7_xy = &dtwnl[is][ngwl*(3+6*ipr7)];
            double *dt7_yz = &dtwnl[is][ngwl*(4+6*ipr7)];
            double *dt7_xz = &dtwnl[is][ngwl*(5+6*ipr7)];
670

671 672 673 674 675 676 677
            double *dt8_xx = &dtwnl[is][ngwl*(0+6*ipr8)];
            double *dt8_yy = &dtwnl[is][ngwl*(1+6*ipr8)];
            double *dt8_zz = &dtwnl[is][ngwl*(2+6*ipr8)];
            double *dt8_xy = &dtwnl[is][ngwl*(3+6*ipr8)];
            double *dt8_yz = &dtwnl[is][ngwl*(4+6*ipr8)];
            double *dt8_xz = &dtwnl[is][ngwl*(5+6*ipr8)];

Francois Gygi committed
678 679 680
            for ( int ig = 0; ig < ngwl; ig++ )
            {
              double v,dv;
Francois Gygi committed
681
              const double tg = kpg[ig];
682

683
              s->dvnlg(l,tg,v,dv);
684

Francois Gygi committed
685 686 687
              const double tgx = kpg_x[ig];
              const double tgy = kpg_y[ig];
              const double tgz = kpg_z[ig];
Francois Gygi committed
688 689 690
              const double tgx2 = tgx * tgx;
              const double tgy2 = tgy * tgy;
              const double tgz2 = tgz * tgz;
691

Francois Gygi committed
692
              const double tgi = kpgi[ig];
693
              const double tg2 = tg * tg;
Francois Gygi committed
694
              const double tgi2 = tgi * tgi;
695

696 697 698 699 700 701
              const double tgxx = tgx2 * tgi2;
              const double tgyy = tgy2 * tgi2;
              const double tgzz = tgz2 * tgi2;
              const double tgxy = tgx * tgy * tgi2;
              const double tgyz = tgy * tgz * tgi2;
              const double tgxz = tgx * tgz * tgi2;
702

703 704 705 706 707
              const double y4 = s54pi * 0.5 * (3.0 * tgzz - 1.0 );
              const double y5 = s54pi * 0.5 * s3 * ( tgxx - tgyy );
              const double y6 = s54pi * s3 * tgxy;
              const double y7 = s54pi * s3 * tgyz;
              const double y8 = s54pi * s3 * tgxz;
708

709 710 711
              const double y1x = s34pi * tgx * tgi;
              const double y1y = s34pi * tgy * tgi;
              const double y1z = s34pi * tgz * tgi;
712

713 714 715 716 717 718
              const double dx_xx = y1x * tgxx - y1x;
              const double dx_yy = y1x * tgyy;
              const double dx_zz = y1x * tgzz;
              const double dx_xy = y1x * tgxy;
              const double dx_yz = y1x * tgyz;
              const double dx_xz = y1x * tgxz;
719

720 721 722 723 724 725
              const double dy_xx = y1y * tgxx;
              const double dy_yy = y1y * tgyy - y1y;
              const double dy_zz = y1y * tgzz;
              const double dy_xy = y1y * tgxy - y1x;
              const double dy_yz = y1y * tgyz;
              const double dy_xz = y1y * tgxz;
726

727 728 729 730 731 732
              const double dz_xx = y1z * tgxx;
              const double dz_yy = y1z * tgyy;
              const double dz_zz = y1z * tgzz - y1z;
              const double dz_xy = y1z * tgxy;
              const double dz_yz = y1z * tgyz - y1y;
              const double dz_xz = y1z * tgxz - y1x;
733

734 735 736 737 738 739
              t4[ig]  = y4 * v;
              t5[ig]  = y5 * v;
              t6[ig]  = y6 * v;
              t7[ig]  = y7 * v;
              t8[ig]  = y8 * v;

Francois Gygi committed
740
              // y4 = s54pi 1/2 ( 3 z^2/r^2 - 1 )
741 742 743 744 745 746
              dt4_xx[ig] = -(v * s20pi * dz_xx * y1z - y4 * dv * tg * tgxx);
              dt4_yy[ig] = -(v * s20pi * dz_yy * y1z - y4 * dv * tg * tgyy);
              dt4_zz[ig] = -(v * s20pi * dz_zz * y1z - y4 * dv * tg * tgzz);
              dt4_xy[ig] = -(v * s20pi * dz_xy * y1z - y4 * dv * tg * tgxy);
              dt4_yz[ig] = -(v * s20pi * dz_yz * y1z - y4 * dv * tg * tgyz);
              dt4_xz[ig] = -(v * s20pi * dz_xz * y1z - y4 * dv * tg * tgxz);
747

Francois Gygi committed
748
              // y5 = s54pi sqrt(3)/2 ( x^2 - y^2 ) / r^2
749 750 751 752 753 754 755
              dt5_xx[ig] = -(v * s20pi3 * (y1x * dx_xx - y1y * dy_xx) - y5 * dv * tg * tgxx);
              dt5_yy[ig] = -(v * s20pi3 * (y1x * dx_yy - y1y * dy_yy) - y5 * dv * tg * tgyy);
              dt5_zz[ig] = -(v * s20pi3 * (y1x * dx_zz - y1y * dy_zz) - y5 * dv * tg * tgzz);
              dt5_xy[ig] = -(v * s20pi3 * (y1x * dx_xy - y1y * dy_xy) - y5 * dv * tg * tgxy);
              dt5_yz[ig] = -(v * s20pi3 * (y1x * dx_yz - y1y * dy_yz) - y5 * dv * tg * tgyz);
              dt5_xz[ig] = -(v * s20pi3 * (y1x * dx_xz - y1y * dy_xz) - y5 * dv * tg * tgxz);

Francois Gygi committed
756
              // y6 = s54pi sqrt(3) x y / r^2
757 758 759 760 761 762 763
              dt6_xx[ig] = -(v * s20pi3 * (dx_xx * y1y + y1x * dy_xx) - y6 * dv * tg * tgxx);
              dt6_yy[ig] = -(v * s20pi3 * (dx_yy * y1y + y1x * dy_yy) - y6 * dv * tg * tgyy);
              dt6_zz[ig] = -(v * s20pi3 * (dx_zz * y1y + y1x * dy_zz) - y6 * dv * tg * tgzz);
              dt6_xy[ig] = -(v * s20pi3 * (dx_xy * y1y + y1x * dy_xy) - y6 * dv * tg * tgxy);
              dt6_yz[ig] = -(v * s20pi3 * (dx_yz * y1y + y1x * dy_yz) - y6 * dv * tg * tgyz);
              dt6_xz[ig] = -(v * s20pi3 * (dx_xz * y1y + y1x * dy_xz) - y6 * dv * tg * tgxz);

Francois Gygi committed
764
              // y7 = s54pi sqrt(3) y z / r^2
765 766 767 768 769 770
              dt7_xx[ig] = -(v * s20pi3 * (dy_xx * y1z + y1y * dz_xx) - y7 * dv * tg * tgxx);
              dt7_yy[ig] = -(v * s20pi3 * (dy_yy * y1z + y1y * dz_yy) - y7 * dv * tg * tgyy);
              dt7_zz[ig] = -(v * s20pi3 * (dy_zz * y1z + y1y * dz_zz) - y7 * dv * tg * tgzz);
              dt7_xy[ig] = -(v * s20pi3 * (dy_xy * y1z + y1y * dz_xy) - y7 * dv * tg * tgxy);
              dt7_yz[ig] = -(v * s20pi3 * (dy_yz * y1z + y1y * dz_yz) - y7 * dv * tg * tgyz);
              dt7_xz[ig] = -(v * s20pi3 * (dy_xz * y1z + y1y * dz_xz) - y7 * dv * tg * tgxz);
771

Francois Gygi committed
772
              // y8 = s54pi sqrt(3) z x / r^2
773 774 775 776 777 778
              dt8_xx[ig] = -(v * s20pi3 * (dx_xx * y1z + y1x * dz_xx) - y8 * dv * tg * tgxx);
              dt8_yy[ig] = -(v * s20pi3 * (dx_yy * y1z + y1x * dz_yy) - y8 * dv * tg * tgyy);
              dt8_zz[ig] = -(v * s20pi3 * (dx_zz * y1z + y1x * dz_zz) - y8 * dv * tg * tgzz);
              dt8_xy[ig] = -(v * s20pi3 * (dx_xy * y1z + y1x * dz_xy) - y8 * dv * tg * tgxy);
              dt8_yz[ig] = -(v * s20pi3 * (dx_yz * y1z + y1x * dz_yz) - y8 * dv * tg * tgyz);
              dt8_xz[ig] = -(v * s20pi3 * (dx_xz * y1z + y1x * dz_xz) - y8 * dv * tg * tgxz);
Francois Gygi committed
779 780 781 782 783 784 785 786 787 788 789 790 791 792 793
            }
          }
          else
          {
            // semi-local
            for ( int iquad = 0; iquad < nquad[is]; iquad++ )
            {
              // twnl[is][ipr][ig]
              // ipr = iquad+nquad[is]*ilm
              // index = ig + ngwl*ipr
              const int ipr4 = iquad+nquad[is]*ilm;
              const int ipr5 = iquad+nquad[is]*(ilm+1);
              const int ipr6 = iquad+nquad[is]*(ilm+2);
              const int ipr7 = iquad+nquad[is]*(ilm+3);
              const int ipr8 = iquad+nquad[is]*(ilm+4);
794

Francois Gygi committed
795 796 797 798 799 800
              double *t4 = &twnl[is][ngwl*ipr4];
              double *t5 = &twnl[is][ngwl*ipr5];
              double *t6 = &twnl[is][ngwl*ipr6];
              double *t7 = &twnl[is][ngwl*ipr7];
              double *t8 = &twnl[is][ngwl*ipr8];

801 802 803 804 805 806 807 808
              // dtwnl[is][ipr][ij][ngwl]
              // index = ig + ngwl * ( ij + 6 * ipr )
              double *dt4_xx = &dtwnl[is][ngwl*(0+6*ipr4)];
              double *dt4_yy = &dtwnl[is][ngwl*(1+6*ipr4)];
              double *dt4_zz = &dtwnl[is][ngwl*(2+6*ipr4)];
              double *dt4_xy = &dtwnl[is][ngwl*(3+6*ipr4)];
              double *dt4_yz = &dtwnl[is][ngwl*(4+6*ipr4)];
              double *dt4_xz = &dtwnl[is][ngwl*(5+6*ipr4)];
Francois Gygi committed
809

810 811 812 813 814 815
              double *dt5_xx = &dtwnl[is][ngwl*(0+6*ipr5)];
              double *dt5_yy = &dtwnl[is][ngwl*(1+6*ipr5)];
              double *dt5_zz = &dtwnl[is][ngwl*(2+6*ipr5)];
              double *dt5_xy = &dtwnl[is][ngwl*(3+6*ipr5)];
              double *dt5_yz = &dtwnl[is][ngwl*(4+6*ipr5)];
              double *dt5_xz = &dtwnl[is][ngwl*(5+6*ipr5)];
816

817 818 819 820 821 822
              double *dt6_xx = &dtwnl[is][ngwl*(0+6*ipr6)];
              double *dt6_yy = &dtwnl[is][ngwl*(1+6*ipr6)];
              double *dt6_zz = &dtwnl[is][ngwl*(2+6*ipr6)];
              double *dt6_xy = &dtwnl[is][ngwl*(3+6*ipr6)];
              double *dt6_yz = &dtwnl[is][ngwl*(4+6*ipr6)];
              double *dt6_xz = &dtwnl[is][ngwl*(5+6*ipr6)];
823

824 825 826 827 828 829
              double *dt7_xx = &dtwnl[is][ngwl*(0+6*ipr7)];
              double *dt7_yy = &dtwnl[is][ngwl*(1+6*ipr7)];
              double *dt7_zz = &dtwnl[is][ngwl*(2+6*ipr7)];
              double *dt7_xy = &dtwnl[is][ngwl*(3+6*ipr7)];
              double *dt7_yz = &dtwnl[is][ngwl*(4+6*ipr7)];
              double *dt7_xz = &dtwnl[is][ngwl*(5+6*ipr7)];
830

831 832 833 834 835 836
              double *dt8_xx = &dtwnl[is][ngwl*(0+6*ipr8)];
              double *dt8_yy = &dtwnl[is][ngwl*(1+6*ipr8)];
              double *dt8_zz = &dtwnl[is][ngwl*(2+6*ipr8)];
              double *dt8_xy = &dtwnl[is][ngwl*(3+6*ipr8)];
              double *dt8_yz = &dtwnl[is][ngwl*(4+6*ipr8)];
              double *dt8_xz = &dtwnl[is][ngwl*(5+6*ipr8)];
Francois Gygi committed
837 838 839 840 841 842 843

              const double r = rquad[is][iquad];
              for ( int ig = 0; ig < ngwl; ig++ )
              {
                double v = 0.0, dv = 0.0;
                // j_2(z) = (3/z^3-1/z) sin(z) - 3/z^2 cos(z)
                // j_2(z) = (1/z)*((3/z^2-1)*sin(z) - (3/z) cos(z) )
844

Francois Gygi committed
845
                const double tg = kpg[ig];
846
                const double z = tg * r;
Francois Gygi committed
847 848 849 850 851 852 853
                if ( z != 0.0 )
                {
                  const double zi = 1.0 / z;
                  const double c = cos(z);
                  const double s = sin(z);
                  const double j2 = ((3.0*zi*zi - 1.0) * s - 3.0*zi * c ) * zi;
                  const double z2 = z * z;
854
                  const double dj2 =
Francois Gygi committed
855 856 857 858 859
                    ( (4.0 * z2 - 9.0) * s + z*(9.0-z2) * c ) / (z2*z2) ;
                  // v = 4 pi j2(Gr) r
                  v = fpi * j2 * r;
                  // dv = d/dG v = 4 pi dj2(Gr)/d(Gr) d(Gr)/dG r
                  //    = 4 pi dj2 r^2
860
                  dv = fpi * dj2 * r * r;
Francois Gygi committed
861
                }
862

Francois Gygi committed
863 864 865
                const double tgx = kpg_x[ig];
                const double tgy = kpg_y[ig];
                const double tgz = kpg_z[ig];
Francois Gygi committed
866 867 868
                const double tgx2 = tgx * tgx;
                const double tgy2 = tgy * tgy;
                const double tgz2 = tgz * tgz;
869

Francois Gygi committed
870
                const double tgi = kpgi[ig];
871
                const double tg2 = tg * tg;
Francois Gygi committed
872
                const double tgi2 = tgi * tgi;
873

874 875 876 877 878 879
                const double tgxx = tgx2 * tgi2;
                const double tgyy = tgy2 * tgi2;
                const double tgzz = tgz2 * tgi2;
                const double tgxy = tgx * tgy * tgi2;
                const double tgyz = tgy * tgz * tgi2;
                const double tgxz = tgx * tgz * tgi2;
880

881 882 883 884 885
                const double y4 = s54pi * 0.5 * (3.0 * tgzz - 1.0 );
                const double y5 = s54pi * 0.5 * s3 * ( tgxx - tgyy );
                const double y6 = s54pi * s3 * tgxy;
                const double y7 = s54pi * s3 * tgyz;
                const double y8 = s54pi * s3 * tgxz;
886

887 888 889
                const double y1x = s34pi * tgx * tgi;
                const double y1y = s34pi * tgy * tgi;
                const double y1z = s34pi * tgz * tgi;
890

891 892 893 894 895 896
                const double dx_xx = y1x * tgxx - y1x;
                const double dx_yy = y1x * tgyy;
                const double dx_zz = y1x * tgzz;
                const double dx_xy = y1x * tgxy;
                const double dx_yz = y1x * tgyz;
                const double dx_xz = y1x * tgxz;
897

898 899 900 901 902 903
                const double dy_xx = y1y * tgxx;
                const double dy_yy = y1y * tgyy - y1y;
                const double dy_zz = y1y * tgzz;
                const double dy_xy = y1y * tgxy - y1x;
                const double dy_yz = y1y * tgyz;
                const double dy_xz = y1y * tgxz;
904

905 906 907 908 909 910
                const double dz_xx = y1z * tgxx;
                const double dz_yy = y1z * tgyy;
                const double dz_zz = y1z * tgzz - y1z;
                const double dz_xy = y1z * tgxy;
                const double dz_yz = y1z * tgyz - y1y;
                const double dz_xz = y1z * tgxz - y1x;
911

Francois Gygi committed
912 913 914 915 916 917 918
                t4[ig]  = y4 * v;
                t5[ig]  = y5 * v;
                t6[ig]  = y6 * v;
                t7[ig]  = y7 * v;
                t8[ig]  = y8 * v;

                // y4 = s54pi 1/2 ( 3 z^2/r^2 - 1 )
919 920 921 922 923 924
                dt4_xx[ig] = -(v * s20pi * dz_xx * y1z - y4 * dv * tg * tgxx);
                dt4_yy[ig] = -(v * s20pi * dz_yy * y1z - y4 * dv * tg * tgyy);
                dt4_zz[ig] = -(v * s20pi * dz_zz * y1z - y4 * dv * tg * tgzz);
                dt4_xy[ig] = -(v * s20pi * dz_xy * y1z - y4 * dv * tg * tgxy);
                dt4_yz[ig] = -(v * s20pi * dz_yz * y1z - y4 * dv * tg * tgyz);
                dt4_xz[ig] = -(v * s20pi * dz_xz * y1z - y4 * dv * tg * tgxz);
925

Francois Gygi committed
926
                // y5 = s54pi sqrt(3)/2 ( x^2 - y^2 ) / r^2
927 928 929 930 931 932 933
                dt5_xx[ig] = -(v * s20pi3 * (y1x * dx_xx - y1y * dy_xx) - y5 * dv * tg * tgxx);
                dt5_yy[ig] = -(v * s20pi3 * (y1x * dx_yy - y1y * dy_yy) - y5 * dv * tg * tgyy);
                dt5_zz[ig] = -(v * s20pi3 * (y1x * dx_zz - y1y * dy_zz) - y5 * dv * tg * tgzz);
                dt5_xy[ig] = -(v * s20pi3 * (y1x * dx_xy - y1y * dy_xy) - y5 * dv * tg * tgxy);
                dt5_yz[ig] = -(v * s20pi3 * (y1x * dx_yz - y1y * dy_yz) - y5 * dv * tg * tgyz);
                dt5_xz[ig] = -(v * s20pi3 * (y1x * dx_xz - y1y * dy_xz) - y5 * dv * tg * tgxz);

Francois Gygi committed
934
                // y6 = s54pi sqrt(3) x y / r^2
935 936 937 938 939 940 941
                dt6_xx[ig] = -(v * s20pi3 * (dx_xx * y1y + y1x * dy_xx) - y6 * dv * tg * tgxx);
                dt6_yy[ig] = -(v * s20pi3 * (dx_yy * y1y + y1x * dy_yy) - y6 * dv * tg * tgyy);
                dt6_zz[ig] = -(v * s20pi3 * (dx_zz * y1y + y1x * dy_zz) - y6 * dv * tg * tgzz);
                dt6_xy[ig] = -(v * s20pi3 * (dx_xy * y1y + y1x * dy_xy) - y6 * dv * tg * tgxy);
                dt6_yz[ig] = -(v * s20pi3 * (dx_yz * y1y + y1x * dy_yz) - y6 * dv * tg * tgyz);
                dt6_xz[ig] = -(v * s20pi3 * (dx_xz * y1y + y1x * dy_xz) - y6 * dv * tg * tgxz);

Francois Gygi committed
942
                // y7 = s54pi sqrt(3) y z / r^2
943 944 945 946 947 948
                dt7_xx[ig] = -(v * s20pi3 * (dy_xx * y1z + y1y * dz_xx) - y7 * dv * tg * tgxx);
                dt7_yy[ig] = -(v * s20pi3 * (dy_yy * y1z + y1y * dz_yy) - y7 * dv * tg * tgyy);
                dt7_zz[ig] = -(v * s20pi3 * (dy_zz * y1z + y1y * dz_zz) - y7 * dv * tg * tgzz);
                dt7_xy[ig] = -(v * s20pi3 * (dy_xy * y1z + y1y * dz_xy) - y7 * dv * tg * tgxy);
                dt7_yz[ig] = -(v * s20pi3 * (dy_yz * y1z + y1y * dz_yz) - y7 * dv * tg * tgyz);
                dt7_xz[ig] = -(v * s20pi3 * (dy_xz * y1z + y1y * dz_xz) - y7 * dv * tg * tgxz);
949

Francois Gygi committed
950
                // y8 = s54pi sqrt(3) z x / r^2
951 952 953 954 955 956
                dt8_xx[ig] = -(v * s20pi3 * (dx_xx * y1z + y1x * dz_xx) - y8 * dv * tg * tgxx);
                dt8_yy[ig] = -(v * s20pi3 * (dx_yy * y1z + y1x * dz_yy) - y8 * dv * tg * tgyy);
                dt8_zz[ig] = -(v * s20pi3 * (dx_zz * y1z + y1x * dz_zz) - y8 * dv * tg * tgzz);
                dt8_xy[ig] = -(v * s20pi3 * (dx_xy * y1z + y1x * dz_xy) - y8 * dv * tg * tgxy);
                dt8_yz[ig] = -(v * s20pi3 * (dx_yz * y1z + y1x * dz_yz) - y8 * dv * tg * tgyz);
                dt8_xz[ig] = -(v * s20pi3 * (dx_xz * y1z + y1x * dz_xz) - y8 * dv * tg * tgxz);
Francois Gygi committed
957 958 959 960 961 962 963 964 965 966 967 968 969
              } // ig
            } // iquad
          }
          ilm += 2*l+1;
        }
        else
        {
          assert(false);
        }
      } // l != lloc[is]
    } // for l
    assert(ilm == s->nlm());
  }
970
  tmap["update_twnl"].stop();
Francois Gygi committed
971 972 973
}

////////////////////////////////////////////////////////////////////////////////
974
double NonLocalPotential::energy(bool compute_hpsi, SlaterDet& dsd,
975
    bool compute_forces, vector<vector<double> >& fion_enl,
976
    bool compute_stress, valarray<double>& sigma_enl)
Francois Gygi committed
977
{
978
  const vector<double>& occ = sd_.occ();
979 980 981 982 983 984 985 986
  const int ngwl = basis_.localsize();
  // define atom block size
  const int na_block_size = 32;
  valarray<double> gr(na_block_size*ngwl); // gr[ig+ia*ngwl]
  valarray<double> cgr(na_block_size*ngwl); // cgr[ig+ia*ngwl]
  valarray<double> sgr(na_block_size*ngwl); // sgr[ig+ia*ngwl]
  vector<vector<double> > tau;
  atoms_.get_positions(tau);
Francois Gygi committed
987 988

  double enl = 0.0;
989
  double tsum[6] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 };
990 991 992
  for ( int is = 0; is < fion_enl.size(); is++ )
    for ( int i = 0; i < fion_enl[is].size(); i++ )
       fion_enl[is][i] = 0.0;
Francois Gygi committed
993 994 995 996 997

  if ( nspnl == 0 ) return 0.0;
  const double omega = basis_.cell().volume();
  assert(omega != 0.0);
  const double omega_inv = 1.0 / omega;
998 999 1000

  for ( int is = 0; is < nsp; is++ )
  {
Francois Gygi committed
1001 1002
    if ( npr[is] > 0 ) // species is is non-local
    {
1003 1004 1005 1006 1007 1008 1009 1010 1011
      valarray<double> tmpfion(3*na[is]);
      tmpfion = 0.0;
      // define number of atom blocks
      const int na_blocks = na[is] / na_block_size +
                            ( na[is] % na_block_size == 0 ? 0 : 1 );

      valarray<double> anl_loc(npr[is]*na_block_size*2*ngwl);
      const int nstloc = sd_.nstloc();
      // fnl_loc[ipra][n]
Francois Gygi committed
1012 1013 1014 1015 1016
      // fnl is real if basis is real, complex otherwise
      const int fnl_loc_size = basis_.real() ? npr[is]*na_block_size*nstloc :
        2*npr[is]*na_block_size*nstloc;
      valarray<double> fnl_loc(fnl_loc_size);
      valarray<double> fnl_buf(fnl_loc_size);
1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033
      for ( int ia_block = 0; ia_block < na_blocks; ia_block++ )
      {
        // process projectors of atoms in block ia_block

        const int iastart = ia_block * na_block_size;
        const int iaend = (ia_block+1) * na_block_size < na[is] ?
                        (ia_block+1) * na_block_size :
                        na[is];
        const int ia_block_size = iaend - iastart;

        // compute cgr[is][ia][ig], sgr[is][ia][ig]
        tmap["comp_eigr"].start();
        int k = 3;
        double mone = -1.0, zero = 0.0;
        char cn='n';

        // next line: const cast is ok since dgemm_ does not modify argument
Francois Gygi committed
1034
        double* kpgx = const_cast<double*>(basis_.kpgx_ptr(0));
1035
        dgemm(&cn,&cn,(int*)&ngwl,(int*)&ia_block_size,&k,&mone,
Francois Gygi committed
1036
              kpgx,(int*)&ngwl, &tau[is][3*iastart],&k,
1037 1038 1039
              &zero,&gr[0],(int*)&ngwl);

        int len = ia_block_size * ngwl;
1040
#if USE_MASSV
1041
        vsincos(&sgr[0],&cgr[0],&gr[0],&len);
1042
#else
1043 1044 1045 1046 1047 1048
        for ( int i = 0; i < len; i++ )
        {
          const double arg = gr[i];
          sgr[i] = sin(arg);
          cgr[i] = cos(arg);
        }
1049
#endif
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 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101
        tmap["comp_eigr"].stop();

        // compute anl_loc
        tmap["comp_anl"].start();
        for ( int ipr = 0; ipr < npr[is]; ipr++ )
        {
          // twnl[is][ig+ngwl*ipr]
          const double * t = &twnl[is][ngwl*ipr];
          const int l = lproj[is][ipr];

          // anl_loc[ig+ipra*ngwl]

          for ( int ia = 0; ia < ia_block_size; ia++ )
          {
            double* a = &anl_loc[2*(ia+ipr*ia_block_size)*ngwl];
            const double* c = &cgr[ia*ngwl];
            const double* s = &sgr[ia*ngwl];
            if ( l == 0 )
            {
              for ( int ig = 0; ig < ngwl; ig++ )
              {
                a[2*ig]   = t[ig] * c[ig];
                a[2*ig+1] = t[ig] * s[ig];
              }
            }
            else if ( l == 1 )
            {
              for ( int ig = 0; ig < ngwl; ig++ )
              {
                /* Next line: -i * eigr */
                /* -i * (a+i*b) = b - i*a */
                a[2*ig]   =  t[ig] * s[ig];
                a[2*ig+1] = -t[ig] * c[ig];
              }
            }
            else if ( l == 2 )
            {
              for ( int ig = 0; ig < ngwl; ig++ )
              {
                // Next line: (-) sign for -eigr
                a[2*ig]   = -t[ig] * c[ig];
                a[2*ig+1] = -t[ig] * s[ig];
              }
            }
          }
        } // ipr
        tmap["comp_anl"].stop();

        // array anl_loc is complete

        // compute fnl[npra][nstloc] = anl^T * c
        int nprnaloc = ia_block_size * npr[is];
Francois Gygi committed
1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128
        if ( basis_.real() )
        {
          double one=1.0;
          char ct='t';
          int twongwl = 2 * ngwl;
          int c_lda = 2*sd_.c().mloc();
          const complex<double>* c = sd_.c().cvalptr();
          tmap["fnl_gemm"].start();
          dgemm(&ct,&cn,&nprnaloc,(int*)&nstloc,&twongwl,&one,
                &anl_loc[0],&twongwl, (double*)c, &c_lda,
                &zero,&fnl_loc[0],&nprnaloc);
          tmap["fnl_gemm"].stop();

          // correct for double counting if ctxt_.myrow() == 0
          if ( ctxt_.myrow() == 0 )
          {
            // rank-one update
            // dger(m,n,alpha,x,incx,y,incy,a,lda);
            // a += alpha * x * transpose(y)
            // x = first row of anl_loc
            // y^T = first row of c
            double alpha = -0.5;
            dger(&nprnaloc,(int*)&nstloc,&alpha,&anl_loc[0],&twongwl,
                 (double*)c,&c_lda,&fnl_loc[0],&nprnaloc);
          }
        }
        else
1129
        {
Francois Gygi committed
1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141
          // fnl is complex
          complex<double> cone=1.0;
          complex<double> czero=0.0;
          char cc='c';
          int c_lda = sd_.c().mloc();
          const complex<double>* c = sd_.c().cvalptr();
          tmap["fnl_zemm"].start();
          zgemm(&cc,&cn,&nprnaloc,(int*)&nstloc,(int*)&ngwl,&cone,
                (complex<double>*) &anl_loc[0],(int*)&ngwl,
                (complex<double>*) c, &c_lda,&czero,
                (complex<double>*)&fnl_loc[0],&nprnaloc);
          tmap["fnl_zemm"].stop();
1142 1143
        }

1144
#if USE_MPI
1145 1146 1147
        tmap["fnl_allreduce"].start();
        // Allreduce fnl partial sum
        MPI_Comm basis_comm = basis_.context().comm();
Francois Gygi committed
1148
        int fnl_size = basis_.real() ? nprnaloc*nstloc : 2*nprnaloc*nstloc;
1149 1150 1151 1152 1153
        MPI_Allreduce(&fnl_loc[0],&fnl_buf[0],fnl_size,
                      MPI_DOUBLE,MPI_SUM,basis_comm);
        tmap["fnl_allreduce"].stop();

        // factor 2.0 in next line is: counting G, -G
Francois Gygi committed
1154 1155 1156 1157
        if ( basis_.real() )
          fnl_loc = 2.0 * fnl_buf;
        else
          fnl_loc = fnl_buf;
1158 1159 1160 1161 1162
#else
        // factor 2.0 in next line is: counting G, -G
        if ( basis_.real() )
          fnl_loc *= 2.0;
#endif
1163 1164 1165

        // accumulate Enl contribution
        const int nbase = ctxt_.mycol() * sd_.c().nb();
Francois Gygi committed
1166
        if ( basis_.real() )
1167
        {
Francois Gygi committed
1168
          for ( int ipr = 0; ipr < npr[is]; ipr++ )
1169
          {
Francois Gygi committed
1170 1171
            const double fac = wt[is][ipr] * omega_inv;
            for ( int n = 0; n < nstloc; n++ )
1172
            {
Francois Gygi committed
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
              const double facn = fac * occ[n + nbase];
              for ( int ia = 0; ia < ia_block_size; ia++ )
              {
                const int i = ia + ipr*ia_block_size + n * nprnaloc;
                //cout << "fnl_loc[ipr=" << ipr << ",ia=" << ia
                //     << ",n=" << n << "]: " << fnl_loc[i] << endl;
                const double tmp = fnl_loc[i];
                enl += facn * tmp * tmp;
                fnl_loc[i] = fac * tmp;
              }
            }
          }
        }
        else
        {
          // fnl is complex
          for ( int ipr = 0; ipr < npr[is]; ipr++ )
          {
            const double fac = wt[is][ipr] * omega_inv;
            for ( int n = 0; n < nstloc; n++ )
            {
              const double facn = fac * occ[n + nbase];
              for ( int ia = 0; ia < ia_block_size; ia++ )
              {
                const int i = ia + ipr*ia_block_size + n * nprnaloc;
                //cout << "fnl_loc[ipr=" << ipr << ",ia=" << ia
                //     << ",n=" << n << "]: " << fnl_loc[i] << endl;
                const double f_re = fnl_loc[2*i];
                const double f_im = fnl_loc[2*i+1];
                enl += facn * (f_re*f_re + f_im*f_im);
                fnl_loc[2*i] = fac * f_re;
                fnl_loc[2*i+1] = fac * f_im;
              }
1206 1207 1208 1209 1210 1211
            }
          }
        }

        if ( compute_hpsi )
        {
Francois Gygi committed
1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237
          if ( basis_.real() )
          {
            tmap["enl_hpsi"].start();
            // compute cp += anl * fnl
            complex<double>* cp = dsd.c().valptr();
            int cp_lda = 2*dsd.c().mloc();
            int twongwl = 2 * ngwl;
            double one = 1.0;
            dgemm(&cn,&cn,&twongwl,(int*)&nstloc,&nprnaloc,&one,
                  &anl_loc[0],&twongwl, &fnl_loc[0],&nprnaloc,
                  &one,(double*)cp, &cp_lda);
            tmap["enl_hpsi"].stop();
          }
          else
          {
            tmap["enl_hpsi"].start();
            // compute cp += anl * fnl
            complex<double>* cp = dsd.c().valptr();
            complex<double> cone = 1.0;
            int cp_lda = dsd.c().mloc();
            zgemm(&cn,&cn,(int*)&ngwl,(int*)&nstloc,&nprnaloc,&cone,
                  (complex<double>*) &anl_loc[0],(int*)&ngwl,
                  (complex<double>*) &fnl_loc[0],&nprnaloc,
                  &cone,cp, &cp_lda);
            tmap["enl_hpsi"].stop();
          }
1238 1239 1240 1241 1242 1243 1244
        }

        // ionic forces
        if ( compute_forces )
        {
          tmap["enl_fion"].start();

Francois Gygi committed
1245
          valarray<double> dfnl_loc(fnl_loc_size);
1246 1247
          for ( int j = 0; j < 3; j++ )
          {
Francois Gygi committed
1248
            const double *const kpgxj = basis_.kpgx_ptr(j);
1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267

            // compute anl_loc
            for ( int ipr = 0; ipr < npr[is]; ipr++ )
            {
              // twnl[is][ig+ngwl*ipr]
              const double * t = &twnl[is][ngwl*ipr];
              const int l = lproj[is][ipr];

              // anl_loc[ig+ipra*ngwl]

              for ( int ia = 0; ia < ia_block_size; ia++ )
              {
                double* a = &anl_loc[2*(ia+ipr*ia_block_size)*ngwl];
                const double* c = &cgr[ia*ngwl];
                const double* s = &sgr[ia*ngwl];
                if ( l == 0 )
                {
                  for ( int ig = 0; ig < ngwl; ig++ )
                  {
Francois Gygi committed
1268
                    const double tt = kpgxj[ig] * t[ig];
1269 1270 1271 1272 1273 1274 1275 1276 1277 1278
                    // Next lines: -i * ( a + ib ) = b - ia
                    a[2*ig]   =  tt * s[ig];
                    a[2*ig+1] = -tt * c[ig];
                  }
                }
                else if ( l == 1 )
                {
                  for ( int ig = 0; ig < ngwl; ig++ )
                  {
                    // Next lines: (-i)**2 * ( a + ib ) = - a - ib
Francois Gygi committed
1279
                    const double tt = - kpgxj[ig] * t[ig];
1280 1281 1282 1283 1284 1285 1286 1287 1288
                    a[2*ig]   = tt * c[ig];
                    a[2*ig+1] = tt * s[ig];
                  }
                }
                else if ( l == 2 )
                {
                  for ( int ig = 0; ig < ngwl; ig++ )
                  {
                    // Next lines: (-i) * - ( a + ib ) = i*(a+ib) = - b + ia
Francois Gygi committed
1289
                    const double tt = kpgxj[ig] * t[ig];
1290 1291 1292 1293 1294 1295 1296 1297 1298 1299
                    a[2*ig]   = -tt * s[ig];
                    a[2*ig+1] =  tt * c[ig];
                  }
                }
              }
            } // ipr

            // array anl_loc is complete

            // compute dfnl[npra][nstloc] = anl^T * c
Francois Gygi committed
1300 1301 1302 1303 1304
            if ( basis_.real() )
            {
              // factor 2.0 in next line is: counting G, -G
              // Note: no need to correct for double counting of the
              // G=0 component which is always zero
1305
              double two=2.0;
Francois Gygi committed
1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327
              char ct='t';
              const int twongwl = 2 * ngwl;
              const int nprnaloc = ia_block_size * npr[is];
              int c_lda = 2*sd_.c().mloc();
              const complex<double>* c = sd_.c().cvalptr();
              dgemm(&ct,&cn,(int*)&nprnaloc,(int*)&nstloc,(int*)&twongwl,&two,
                    &anl_loc[0],(int*)&twongwl, (double*)c,(int*)&c_lda,
                    &zero,&dfnl_loc[0],(int*)&nprnaloc);
            }
            else
            {
              complex<double> cone=1.0;
              complex<double> czero=0.0;
              char cc='c';
              const int nprnaloc = ia_block_size * npr[is];
              int c_lda = sd_.c().mloc();
              const complex<double>* c = sd_.c().cvalptr();
              zgemm(&cc,&cn,(int*)&nprnaloc,(int*)&nstloc,(int*)&ngwl,&cone,
                    (complex<double>*) &anl_loc[0],(int*)&ngwl,
                    (complex<double>*) c,(int*)&c_lda,
                    &czero,(complex<double>*)&dfnl_loc[0],(int*)&nprnaloc);
            }
1328 1329

            // accumulate non-local contributions to forces
Francois Gygi committed
1330
            if ( basis_.real() )
1331
            {
Francois Gygi committed
1332
              for ( int ipr = 0; ipr < npr[is]; ipr++ )
1333
              {
Francois Gygi committed
1334
                for ( int n = 0; n < nstloc; n++ )
1335
                {
Francois Gygi committed
1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370
                  // Factor 2.0 in next line from derivative of |Fnl|^2
                  const double facn = 2.0 * occ[n + nbase];
                  for ( int ia = 0; ia < ia_block_size; ia++ )
                  {
                    const int ia_global = ia + iastart;
                    const int i = ia + ipr*ia_block_size + n * nprnaloc;
                    //cout << "fnl_loc[ipr=" << ipr << ",ia=" << ia
                    //     << ",n=" << n << "]: " << fnl_loc[i] << endl;
                    tmpfion[3*ia_global+j] -= facn *
                                              fnl_loc[i] * dfnl_loc[i];
                  }
                }
              }
            }
            else
            {
              for ( int ipr = 0; ipr < npr[is]; ipr++ )
              {
                for ( int n = 0; n < nstloc; n++ )
                {
                  // Factor 2.0 in next line from derivative of |Fnl|^2
                  const double facn = 2.0 * occ[n + nbase];
                  for ( int ia = 0; ia < ia_block_size; ia++ )
                  {
                    const int ia_global = ia + iastart;
                    const int i = ia + ipr*ia_block_size + n * nprnaloc;
                    //cout << "fnl_loc[ipr=" << ipr << ",ia=" << ia
                    //     << ",n=" << n << "]: " << fnl_loc[i] << endl;
                    double f_re = fnl_loc[2*i];
                    double f_im = fnl_loc[2*i+1];
                    double df_re = dfnl_loc[2*i];
                    double df_im = dfnl_loc[2*i+1];
                    tmpfion[3*ia_global+j] -=
                      facn * ( f_re * df_re + f_im * df_im );
                  }
1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381
                }
              }
            }
          } // j

          tmap["enl_fion"].stop();
        } // compute_forces

        if ( compute_stress )
        {
          tmap["enl_sigma"].start();
Francois Gygi committed
1382
          valarray<double> dfnl_loc(fnl_loc_size);
1383 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 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455

          for ( int ij = 0; ij < 6; ij++ )
          {
            // compute anl_loc
            int ipr = 0;
            while ( ipr < npr[is] )
            {
              // twnl[is][ig+ngwl*ipr]
              const int l = lproj[is][ipr];
              if ( l == 0 )
              {
                // dtwnl[is][ipr][ij][ngwl]
                // index = ig + ngwl * ( ij + 6 * ipr))
                // ipr = iquad + nquad[is] * ilm, where ilm = 0
                const double *const dt0 = &dtwnl[is][ngwl*(ij+6*ipr)];
                for ( int ia = 0; ia < ia_block_size; ia++ )
                {
                  double* a0 = &anl_loc[2*(ia+ipr*ia_block_size)*ngwl];
                  const double* c = &cgr[ia*ngwl];
                  const double* s = &sgr[ia*ngwl];
                  for ( int ig = 0; ig < ngwl; ig++ )
                  {
                    const double d0 = dt0[ig];
                    // Next lines: -i * ( a + ib ) = b - ia
                    a0[2*ig]   =  d0 * c[ig];
                    a0[2*ig+1] =  d0 * s[ig];
                  }
                }
              }
              else if ( l == 1 )
              {
                const int ipr1 = ipr;
                const int ipr2 = ipr + 1;
                const int ipr3 = ipr + 2;
                // dtwnl[is][ipr][ij][ngwl]
                // index = ig + ngwl * ( ij + 6 * iprx ))
                const double *dt1 = &dtwnl[is][ngwl*(ij+6*ipr1)];
                const double *dt2 = &dtwnl[is][ngwl*(ij+6*ipr2)];
                const double *dt3 = &dtwnl[is][ngwl*(ij+6*ipr3)];
                for ( int ia = 0; ia < ia_block_size; ia++ )
                {
                  double* a1 = &anl_loc[2*(ia+ipr1*ia_block_size)*ngwl];
                  double* a2 = &anl_loc[2*(ia+ipr2*ia_block_size)*ngwl];
                  double* a3 = &anl_loc[2*(ia+ipr3*ia_block_size)*ngwl];
                  const double* c = &cgr[ia*ngwl];
                  const double* s = &sgr[ia*ngwl];
                  for ( int ig = 0; ig < ngwl; ig++ )
                  {
                    const double d1 = dt1[ig];
                    const double d2 = dt2[ig];
                    const double d3 = dt3[ig];
                    // Next line: (-i)^l factor is -i
                    // Next line: -i * eigr
                    // -i * (a+i*b) = b - i*a
                    const double tc = -c[ig]; //  -cosgr[ia][ig]
                    const double ts =  s[ig]; //   singr[ia][ig]
                    a1[2*ig]   = d1 * ts;
                    a1[2*ig+1] = d1 * tc;
                    a2[2*ig]   = d2 * ts;
                    a2[2*ig+1] = d2 * tc;
                    a3[2*ig]   = d3 * ts;
                    a3[2*ig+1] = d3 * tc;
                  }
                }
              }
              else if ( l == 2 )
              {
                const int ipr4 = ipr;
                const int ipr5 = ipr + 1;
                const int ipr6 = ipr + 2;
                const int ipr7 = ipr + 3;
                const int ipr8 = ipr + 4;
                // dtwnl[is][ipr][iquad][ij][ngwl]
1456
                // index = ig + ngwl * ( ij + 6 * ( iquad + nquad[is] * ipr ))
1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504
                const double *dt4 = &dtwnl[is][ngwl*(ij+6*ipr4)];
                const double *dt5 = &dtwnl[is][ngwl*(ij+6*ipr5)];
                const double *dt6 = &dtwnl[is][ngwl*(ij+6*ipr6)];
                const double *dt7 = &dtwnl[is][ngwl*(ij+6*ipr7)];
                const double *dt8 = &dtwnl[is][ngwl*(ij+6*ipr8)];
                for ( int ia = 0; ia < ia_block_size; ia++ )
                {
                  double* a4 = &anl_loc[2*(ia+ipr4*ia_block_size)*ngwl];
                  double* a5 = &anl_loc[2*(ia+ipr5*ia_block_size)*ngwl];
                  double* a6 = &anl_loc[2*(ia+ipr6*ia_block_size)*ngwl];
                  double* a7 = &anl_loc[2*(ia+ipr7*ia_block_size)*ngwl];
                  double* a8 = &anl_loc[2*(ia+ipr8*ia_block_size)*ngwl];
                  const double* c = &cgr[ia*ngwl];
                  const double* s = &sgr[ia*ngwl];

                  for ( int ig = 0; ig < ngwl; ig++ )
                  {
                    const double d4 = dt4[ig];
                    const double d5 = dt5[ig];
                    const double d6 = dt6[ig];
                    const double d7 = dt7[ig];
                    const double d8 = dt8[ig];
                    // Next lines: (-i)^2 * ( a + ib ) =  - ( a + ib )
                    const double tc = -c[ig]; //  -cosgr[ia][ig]
                    const double ts = -s[ig]; //  -singr[ia][ig]
                    a4[2*ig]   = d4 * tc;
                    a4[2*ig+1] = d4 * ts;
                    a5[2*ig]   = d5 * tc;
                    a5[2*ig+1] = d5 * ts;
                    a6[2*ig]   = d6 * tc;
                    a6[2*ig+1] = d6 * ts;
                    a7[2*ig]   = d7 * tc;
                    a7[2*ig+1] = d7 * ts;