NonLocalPotential.C 81.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 18 19
// NonLocalPotential.C
//
////////////////////////////////////////////////////////////////////////////////

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

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


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

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

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

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

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

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

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

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

Francois Gygi committed
90 91 92 93 94 95 96 97 98 99 100 101
      // 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];
102

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

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

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

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

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
      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
161
        //cout << " NonLocalPotential::init: trapezoidal rule, right endpoint"
162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181
        //     << 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
182
      {
183
        assert(false);
Francois Gygi committed
184
      }
185

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

      // compute lproj[is][ipr]
Francois Gygi committed
190
      int ipr_base = 0;
191
      for ( int iop = 0; iop < nop[is]; iop++ )
Francois Gygi committed
192
      {
193
        int l = s->l(iop);
Francois Gygi committed
194 195 196 197 198 199 200 201 202 203
        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;
204
              wt[is][ipr] = s->wsg(iop);
Francois Gygi committed
205 206 207 208 209 210 211 212 213 214
              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

244 245
  if ( nspnl == 0 ) return;

Francois Gygi committed
246 247 248 249
  const int ngwl = basis_.localsize();
  const double pi = M_PI;
  const double fpi = 4.0 * pi;
  const double s14pi = sqrt(1.0/fpi);
250
  const double s34pi = sqrt(3.0/fpi);
Francois Gygi committed
251
  const double s54pi = sqrt(5.0/fpi);
252
  const double s74pi = sqrt(7.0/fpi);
253 254
  const double s20pi = sqrt(20.0*pi);
  const double s20pi3 = sqrt(20.0*pi/3.0);
Francois Gygi committed
255
  const double s3 = sqrt(3.0);
256 257
  const double s32 = sqrt(1.5);
  const double s52 = sqrt(2.5);
258
  const double s15 = sqrt(15.0);
Francois Gygi committed
259

Francois Gygi committed
260 261 262 263 264
  const double *kpg   = basis_.kpg_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);
265

Francois Gygi committed
266 267 268 269
  // compute twnl and dtwnl
  for ( int is = 0; is < nsp; is++ )
  {
    Species *s = atoms_.species_list[is];
270

Francois Gygi committed
271
    int ilm = 0;
272
    for ( int iop = 0; iop < nop[is]; iop++ )
273
    {
274
      int l = s->l(iop);
Francois Gygi committed
275 276 277 278 279 280 281
      if ( l != lloc[is] )
      {
        if ( l == 0 )
        {
          if ( nquad[is] == 0 )
          {
            // Kleinman-Bylander
282

Francois Gygi committed
283
            // twnl[is][ipr][ig]
284
            const int ipr = ilm;
Francois Gygi committed
285
            // index = ig + ngwl*ipr, i.e. index = ig
286
            double *t0 = &twnl[is][ngwl*ipr];
287

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

              t0[ig] = s14pi * v;
304

Francois Gygi committed
305 306 307
              const double tgx = kpg_x[ig];
              const double tgy = kpg_y[ig];
              const double tgz = kpg_z[ig];
308

Francois Gygi committed
309
              const double tmp = kpgi[ig] * s14pi * dv;
310 311 312 313 314 315
              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
316 317 318 319 320 321 322 323 324 325 326 327 328
            }
          }
          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]
329
              // index = ig + ngwl * ( ij + 6 * iquad)
330 331 332 333 334 335
              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
336
              const double r = rquad[is][iquad];
Francois Gygi committed
337 338 339

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

350 351
                // next line: sin(Gr)/(Gr) replaced by 1
                t0[ig] = fpi * s14pi * r;
352

353
                // dtwnl = fpi s14pi G_i G_j / G (r cos(Gr)/G -sin(Gr)/G^2)
354 355 356 357 358 359 360
                // 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;
361 362 363
              }
              else
              {
Francois Gygi committed
364
                // Use the normal procedure
365 366 367 368 369
                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
370
                const double arg = kpg[ig] * r;
371
                // Note: for G=0, gi[0] = 0
372

Francois Gygi committed
373 374 375 376
                const double tgx = kpg_x[ig];
                const double tgy = kpg_y[ig];
                const double tgz = kpg_z[ig];
                const double tgi = kpgi[ig];
377
                const double tgi2 = tgi * tgi;
378

379 380
                const double ts = sin(arg);
                const double tc = cos(arg);
381

382
                t0[ig] = fpi * s14pi * ts * tgi;
383

Francois Gygi committed
384 385
                // dtwnl = fpi s14pi k+G_i k+G_j / k+G
                // (r cos(k+Gr)/k+G -sin(k+Gr)/k+G^2)
386 387 388 389 390 391 392 393 394 395
                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
396
              {
Francois Gygi committed
397 398 399
                // 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
400
                // Ylm = s14pi
Francois Gygi committed
401
                const double arg = kpg[ig] * r;
402

Francois Gygi committed
403 404 405 406
                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
407
                const double tgi2 = tgi * tgi;
408

Francois Gygi committed
409 410
                const double ts = sin(arg);
                const double tc = cos(arg);
411

Francois Gygi committed
412
                t0[ig] = fpi * s14pi * ts * tgi;
413

Francois Gygi committed
414 415
                // 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)
416 417 418 419 420 421 422
                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
423
              }
Francois Gygi committed
424
            } // for iquad
Francois Gygi committed
425 426 427 428 429 430 431 432
          }
          ilm += 2*l+1;
        }
        else if ( l == 1 )
        {
          if ( nquad[is] == 0 )
          {
            // Kleinman-Bylander
433

Francois Gygi committed
434 435 436 437 438 439 440 441 442
            // 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];
443

444 445 446 447 448 449 450 451 452 453 454 455 456 457 458
            // 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
459

460 461 462 463 464 465
            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
466 467 468 469

            for ( int ig = 0; ig < ngwl; ig++ )
            {
              double v,dv;
Francois Gygi committed
470
              const double tg = kpg[ig];
471
              s->dvnlg(iop,tg,v,dv);
472

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

Francois Gygi committed
480
              const double tgi = kpgi[ig];
Francois Gygi committed
481
              const double tgi2 = tgi * tgi;
482

Francois Gygi committed
483 484 485
              const double y1 = s34pi * tgx * tgi;
              const double y2 = s34pi * tgy * tgi;
              const double y3 = s34pi * tgz * tgi;
486

Francois Gygi committed
487 488 489 490
              t1[ig]  = y1 * v;
              t2[ig]  = y2 * v;
              t3[ig]  = y3 * v;

491 492 493 494 495 496 497 498
              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;
499

500 501 502 503 504 505 506 507
              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;
508

509 510 511 512 513 514 515 516
              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
517 518 519 520 521 522 523 524 525 526 527 528 529 530 531
            }
          }
          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];
532

Francois Gygi committed
533
              // dtwnl[is][ipr][j][ngwl]
534 535 536 537 538 539 540
              // 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
541

542 543 544 545 546 547
              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
548

549 550 551 552 553 554
              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)];
555

Francois Gygi committed
556 557 558 559 560
              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
561
                const double tg = kpg[ig];
562
                const double z = tg * r;
Francois Gygi committed
563 564 565 566 567 568
                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;
569
                  const double dj1 =
Francois Gygi committed
570 571 572 573 574 575 576
                    ( 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;
                }
577

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

Francois Gygi committed
585
                const double tgi = kpgi[ig];
Francois Gygi committed
586
                const double tgi2 = tgi * tgi;
587

Francois Gygi committed
588 589 590
                const double y1 = s34pi * tgx * tgi;
                const double y2 = s34pi * tgy * tgi;
                const double y3 = s34pi * tgz * tgi;
591

Francois Gygi committed
592 593 594 595
                t1[ig]  = y1 * v;
                t2[ig]  = y2 * v;
                t3[ig]  = y3 * v;

596 597 598 599 600 601 602 603
                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;
604

605 606 607 608 609 610 611 612
                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;
613

614 615 616 617 618 619 620 621
                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
622 623 624 625 626 627 628 629 630 631 632 633 634 635 636
              } // 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;
637

Francois Gygi committed
638 639 640 641 642
            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];
643

644 645 646 647 648 649 650 651
            // 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
652

653 654 655 656 657 658
            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)];
659

660 661 662 663 664 665
            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)];
666

667 668 669 670 671 672
            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)];
673

674 675 676 677 678 679 680
            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
681 682 683
            for ( int ig = 0; ig < ngwl; ig++ )
            {
              double v,dv;
Francois Gygi committed
684
              const double tg = kpg[ig];
685

686
              s->dvnlg(iop,tg,v,dv);
687

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

Francois Gygi committed
695
              const double tgi = kpgi[ig];
Francois Gygi committed
696
              const double tgi2 = tgi * tgi;
697

698 699 700 701 702 703
              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;
704

705 706 707 708 709
              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;
710

711 712 713
              const double y1x = s34pi * tgx * tgi;
              const double y1y = s34pi * tgy * tgi;
              const double y1z = s34pi * tgz * tgi;
714

715 716 717 718 719 720
              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;
721

722 723 724 725 726 727
              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;
728

729 730 731 732 733 734
              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;
735

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

Francois Gygi committed
742
              // y4 = s54pi 1/2 ( 3 z^2/r^2 - 1 )
743 744 745 746 747 748
              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);
749

Francois Gygi committed
750
              // y5 = s54pi sqrt(3)/2 ( x^2 - y^2 ) / r^2
751 752 753 754 755 756 757
              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
758
              // y6 = s54pi sqrt(3) x y / r^2
759 760 761 762 763 764 765
              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
766
              // y7 = s54pi sqrt(3) y z / r^2
767 768 769 770 771 772
              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);
773

Francois Gygi committed
774
              // y8 = s54pi sqrt(3) z x / r^2
775 776 777 778 779 780
              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
781 782 783 784 785 786 787 788 789 790 791 792 793 794 795
            }
          }
          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);
796

Francois Gygi committed
797 798 799 800 801 802
              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];

803 804 805 806 807 808 809 810
              // 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
811

812 813 814 815 816 817
              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)];
818

819 820 821 822 823 824
              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)];
825

826 827 828 829 830 831
              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)];
832

833 834 835 836 837 838
              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
839 840 841 842 843 844 845

              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) )
846

Francois Gygi committed
847
                const double tg = kpg[ig];
848
                const double z = tg * r;
Francois Gygi committed
849 850 851 852 853 854 855
                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;
856
                  const double dj2 =
Francois Gygi committed
857 858 859 860 861
                    ( (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
862
                  dv = fpi * dj2 * r * r;
Francois Gygi committed
863
                }
864

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

Francois Gygi committed
872
                const double tgi = kpgi[ig];
Francois Gygi committed
873
                const double tgi2 = tgi * tgi;
874

875 876 877 878 879 880
                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;
881

882 883 884 885 886
                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;
887

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

892 893 894 895 896 897
                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;
898

899 900 901 902 903 904
                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;
905

906 907 908 909 910 911
                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;
912

Francois Gygi committed
913 914 915 916 917 918 919
                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 )
920 921 922 923 924 925
                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);
926

Francois Gygi committed
927
                // y5 = s54pi sqrt(3)/2 ( x^2 - y^2 ) / r^2
928 929 930 931 932 933 934
                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
935
                // y6 = s54pi sqrt(3) x y / r^2
936 937 938 939 940 941 942
                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
943
                // y7 = s54pi sqrt(3) y z / r^2
944 945 946 947 948 949
                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);
950

Francois Gygi committed
951
                // y8 = s54pi sqrt(3) z x / r^2
952 953 954 955 956 957
                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
958 959 960 961 962
              } // ig
            } // iquad
          }
          ilm += 2*l+1;
        }
963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 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 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 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 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155
        else if ( l == 3 )
        {
          // only implemented for Kleiman-Bylander type
          assert(nquad[is] == 0);
          // Kleinman-Bylander
          const int ipr09 = ilm;
          const int ipr10 = ilm + 1;
          const int ipr11 = ilm + 2;
          const int ipr12 = ilm + 3;
          const int ipr13 = ilm + 4;
          const int ipr14 = ilm + 5;
          const int ipr15 = ilm + 6;

          double *t09 = &twnl[is][ngwl * ipr09];
          double *t10 = &twnl[is][ngwl * ipr10];
          double *t11 = &twnl[is][ngwl * ipr11];
          double *t12 = &twnl[is][ngwl * ipr12];
          double *t13 = &twnl[is][ngwl * ipr13];
          double *t14 = &twnl[is][ngwl * ipr14];
          double *t15 = &twnl[is][ngwl * ipr15];

          // dtwnl[is][ipr][ij][ngwl]
          // index = ig + ngwl * ( ij + 6 * ipr )
          double *dt09_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr09 )];
          double *dt09_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr09 )];
          double *dt09_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr09 )];
          double *dt09_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr09 )];
          double *dt09_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr09 )];
          double *dt09_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr09 )];

          double *dt10_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr10 )];
          double *dt10_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr10 )];
          double *dt10_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr10 )];
          double *dt10_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr10 )];
          double *dt10_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr10 )];
          double *dt10_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr10 )];

          double *dt11_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr11 )];
          double *dt11_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr11 )];
          double *dt11_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr11 )];
          double *dt11_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr11 )];
          double *dt11_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr11 )];
          double *dt11_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr11 )];

          double *dt12_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr12 )];
          double *dt12_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr12 )];
          double *dt12_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr12 )];
          double *dt12_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr12 )];
          double *dt12_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr12 )];
          double *dt12_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr12 )];

          double *dt13_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr13 )];
          double *dt13_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr13 )];
          double *dt13_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr13 )];
          double *dt13_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr13 )];
          double *dt13_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr13 )];
          double *dt13_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr13 )];

          double *dt14_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr14 )];
          double *dt14_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr14 )];
          double *dt14_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr14 )];
          double *dt14_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr14 )];
          double *dt14_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr14 )];
          double *dt14_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr14 )];

          double *dt15_xx = &dtwnl[is][ngwl * ( 0 + 6 * ipr15 )];
          double *dt15_yy = &dtwnl[is][ngwl * ( 1 + 6 * ipr15 )];
          double *dt15_zz = &dtwnl[is][ngwl * ( 2 + 6 * ipr15 )];
          double *dt15_xy = &dtwnl[is][ngwl * ( 3 + 6 * ipr15 )];
          double *dt15_yz = &dtwnl[is][ngwl * ( 4 + 6 * ipr15 )];
          double *dt15_xz = &dtwnl[is][ngwl * ( 5 + 6 * ipr15 )];

          for ( int ig = 0; ig < ngwl; ig++ )
          {
            double v, dv;
            const double tg = kpg[ig];

            s->dvnlg(iop,tg,v,dv);

            const double tgx = kpg_x[ig];
            const double tgy = kpg_y[ig];
            const double tgz = kpg_z[ig];
            const double tgx2 = tgx * tgx;
            const double tgy2 = tgy * tgy;
            const double tgz2 = tgz * tgz;
            const double tgx3 = tgx * tgx2;
            const double tgy3 = tgy * tgy2;
            const double tgz3 = tgz * tgz2;

            const double tgi = kpgi[ig];
            const double tgi2 = tgi * tgi;
            const double tgi3 = tgi * tgi2;

            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;

            const double tgxxx = tgx3 * tgi3;
            const double tgyyy = tgy3 * tgi3;
            const double tgzzz = tgz3 * tgi3;
            const double tgxyy = tgx * tgy2 * tgi3;
            const double tgxzz = tgx * tgz2 * tgi3;
            const double tgyxx = tgy * tgx2 * tgi3;
            const double tgyzz = tgy * tgz2 * tgi3;
            const double tgzxx = tgz * tgx2 * tgi3;
            const double tgzyy = tgz * tgy2 * tgi3;
            const double tgxyz = tgx * tgy * tgz * tgi3;

            const double factor = 0.5 * s74pi;

            const double y09 = factor * s52 * ( 3.0 * tgyxx - tgyyy );
            const double y10 = factor * s15 * 2.0 * tgxyz;
            const double y11 = factor * s32 * ( 4.0 * tgyzz - tgyxx - tgyyy );
            const double y12 = factor
              * ( 2.0 * tgzzz - 3.0 * ( tgzxx + tgzyy ) );
            const double y13 = factor * s32 * ( 4.0 * tgxzz - tgxxx - tgxyy );
            const double y14 = factor * s15 * ( tgzxx - tgzyy );
            const double y15 = factor * s52 * ( tgxxx - 3.0 * tgxyy );

            // derivative of x^3/r^3 w.r.t. x, y, and z
            const double dx_xxx = 3.0 * tgxx * ( 1.0 - tgxx ) * tgi;
            const double dy_xxx = -3.0 * tgxx * tgxy * tgi;
            const double dz_xxx = -3.0 * tgxx * tgxz * tgi;
            // derivative of y^3/r^3 w.r.t. x, y, and z
            const double dx_yyy = -3.0 * tgyy * tgxy * tgi;
            const double dy_yyy = 3.0 * tgyy * ( 1.0 - tgyy ) * tgi;
            const double dz_yyy = -3.0 * tgyy * tgyz * tgi;
            // derivative of z^3/r^3 w.r.t. x, y, and z
            const double dx_zzz = -3.0 * tgzz * tgxz * tgi;
            const double dy_zzz = -3.0 * tgzz * tgyz * tgi;
            const double dz_zzz = 3.0 * tgzz * ( 1.0 - tgzz ) * tgi;
            // derivative of xy^2/r^3 w.r.t. x, y, and z
            const double dx_xyy = tgyy * ( 1.0 - 3.0 * tgxx ) * tgi;
            const double dy_xyy = tgxy * ( 2.0 - 3.0 * tgyy ) * tgi;
            const double dz_xyy = -3.0 * tgyy * tgxz * tgi;
            // derivative of xz^2/r^3 w.r.t. x, y, and z
            const double dx_xzz = tgzz * ( 1.0 - 3.0 * tgxx ) * tgi;
            const double dy_xzz = -3.0 * tgzz * tgxy * tgi;
            const double dz_xzz = tgxz * ( 2.0 - 3.0 * tgzz ) * tgi;
            // derivative of yx^2/r^3 w.r.t. x, y, and z
            const double dx_yxx = tgxy * ( 2.0 - 3.0 * tgxx ) * tgi;
            const double dy_yxx = tgxx * ( 1.0 - 3.0 * tgyy ) * tgi;
            const double dz_yxx = -3.0 * tgxx * tgyz * tgi;
            // derivative of yz^2/r^3 w.r.t. x, y, and z
            const double dx_yzz = -3.0 * tgzz * tgxy * tgi;
            const double dy_yzz = tgzz * ( 1.0 - 3.0 * tgyy ) * tgi;
            const double dz_yzz = tgyz * ( 2.0 - 3.0 * tgzz ) * tgi;
            // derivative of zx^2/r^3 w.r.t. x, y, and z
            const double dx_zxx = tgxz * ( 2.0 - 3.0 * tgxx ) * tgi;
            const double dy_zxx = -3.0 * tgxx * tgyz * tgi;
            const double dz_zxx = tgxx * ( 1.0 - 3.0 * tgzz ) * tgi;
            // derivative of zy^2/r^3 w.r.t. x, y, and z
            const double dx_zyy = -3.0 * tgyy * tgxz * tgi;
            const double dy_zyy = tgyz * ( 2.0 - 3.0 * tgyy ) * tgi;
            const double dz_zyy = tgyy * ( 1.0 - 3.0 * tgzz ) * tgi;
            // derivative of xyz/r^3 w.r.t. x, y, and z
            const double dx_xyz = tgyz * ( 1 - 3.0 * tgxx ) * tgi;
            const double dy_xyz = tgxz * ( 1 - 3.0 * tgyy ) * tgi;
            const double dz_xyz = tgxy * ( 1 - 3.0 * tgzz ) * tgi;

            // derivatives of spherical harmonics

            // y9 = factor * s52 * ( 3.0 * tgyxx - tgyyy );
            const double dx_y09 = factor * s52 * ( 3.0 * dx_yxx - dx_yyy );
            const double dy_y09 = factor * s52 * ( 3.0 * dy_yxx - dy_yyy );
            const double dz_y09 = factor * s52 * ( 3.0 * dz_yxx - dz_yyy );
            // y10 = factor * s15 * 2.0 * tgxyz;
            const double dx_y10 = factor * s15 * 2.0 * dx_xyz;
            const double dy_y10 = factor * s15 * 2.0 * dy_xyz;
            const double dz_y10 = factor * s15 * 2.0 * dz_xyz;
            // y11 = factor * s32 * ( 4.0 * tgyzz - tgyxx - tgyyy );
            const double dx_y11 = factor * s32 * ( 4.0 * dx_yzz - dx_yxx - dx_yyy );
            const double dy_y11 = factor * s32 * ( 4.0 * dy_yzz - dy_yxx - dy_yyy );
            const double dz_y11 = factor * s32 * ( 4.0 * dz_yzz - dz_yxx - dz_yyy );
            // y12 = factor * ( 2.0 * tgzzz - 3.0 * ( tgzxx + tgzyy ) );
            const double dx_y12 = factor
              * ( 2.0 * dx_zzz - 3.0 * ( dx_zxx + dx_zyy ) );
            const double dy_y12 = factor
              * ( 2.0 * dy_zzz - 3.0 * ( dy_zxx + dy_zyy ) );
            const double dz_y12 = factor
              * ( 2.0 * dz_zzz - 3.0 * ( dz_zxx + dz_zyy ) );
            // y13 = factor * s32 * ( 4.0 * tgxzz - tgxxx - tgxyy );
            const double dx_y13 = factor * s32 * ( 4.0 * dx_xzz - dx_xxx - dx_xyy );
            const double dy_y13 = factor * s32 * ( 4.0 * dy_xzz - dy_xxx - dy_xyy );
            const double dz_y13 = factor * s32 * ( 4.0 * dz_xzz - dz_xxx - dz_xyy );
            // y14 = factor * s15 * ( tgzxx - tgzyy );
            const double dx_y14 = factor * s15 * ( dx_zxx - dx_zyy );
            const double dy_y14 = factor * s15 * ( dy_zxx - dy_zyy );
            const double dz_y14 = factor * s15 * ( dz_zxx - dz_zyy );
            // y15 = factor * s52 * ( tgxxx - 3.0 * tgxyy );
1156
            const double dx_y15 = factor * s52 * ( dx_xxx - 3.0 * dx_xyy );
1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173
            const double dy_y15 = factor * s52 * ( dy_xxx - 3.0 * dy_xyy );
            const double dz_y15 = factor * s52 * ( dz_xxx - 3.0 * dz_xyy );

            t09[ig] = y09 * v;
            t10[ig] = y10 * v;
            t11[ig] = y11 * v;
            t12[ig] = y12 * v;
            t13[ig] = y13 * v;
            t14[ig] = y14 * v;
            t15[ig] = y15 * v;

            // contribution to stress tensor
            // sigma_ij = 0.5 * (xi * dj + xj * di)

            dt09_xx[ig] = -( v * dx_y09 * tgx - y09 * dv * tg * tgxx );
            dt09_yy[ig] = -( v * dy_y09 * tgy - y09 * dv * tg * tgyy );
            dt09_zz[ig] = -( v * dz_y09 * tgz - y09 * dv * tg * tgzz );
1174
            dt09_xy[ig] = -( 0.5 * v * ( dx_y09 * tgy + dy_y09 * tgx )
1175
              - y09 * dv * tg * tgxy );
1176
            dt09_yz[ig] = -( 0.5 * v * ( dy_y09 * tgz + dz_y09 * tgy )
1177
              - y09 * dv * tg * tgyz );
1178
            dt09_xz[ig] = -( 0.5 * v * ( dx_y09 * tgz + dz_y09 * tgx )
1179 1180 1181 1182 1183
              - y09 * dv * tg * tgxz );

            dt10_xx[ig] = -( v * dx_y10 * tgx - y10 * dv * tg * tgxx );
            dt10_yy[ig] = -( v * dy_y10 * tgy - y10 * dv * tg * tgyy );
            dt10_zz[ig] = -( v * dz_y10 * tgz - y10 * dv * tg * tgzz );
1184
            dt10_xy[ig] = -( 0.5 * v * ( dx_y10 * tgy + dy_y10 * tgx )
1185
              - y10 * dv * tg * tgxy );
1186
            dt10_yz[ig] = -( 0.5 * v * ( dy_y10 * tgz + dz_y10 * tgy )
1187
              - y10 * dv * tg * tgyz );
1188
            dt10_xz[ig] = -( 0.5 * v * ( dx_y10 * tgz + dz_y10 * tgx )
1189 1190 1191 1192 1193
              - y10 * dv * tg * tgxz );

            dt11_xx[ig] = -( v * dx_y11 * tgx - y11 * dv * tg * tgxx );
            dt11_yy[ig] = -( v * dy_y11 * tgy - y11 * dv * tg * tgyy );
            dt11_zz[ig] = -( v * dz_y11 * tgz - y11 * dv * tg * tgzz );
1194
            dt11_xy[ig] = -( 0.5 * v * ( dx_y11 * tgy + dy_y11 * tgx )
1195
              - y11 * dv * tg * tgxy );
1196
            dt11_yz[ig] = -( 0.5 * v * ( dy_y11 * tgz + dz_y11 * tgy )
1197
              - y11 * dv * tg * tgyz );
1198
            dt11_xz[ig] = -( 0.5 * v * ( dx_y11 * tgz + dz_y11 * tgx )
1199 1200 1201 1202 1203
              - y11 * dv * tg * tgxz );

            dt12_xx[ig] = -( v * dx_y12 * tgx - y12 * dv * tg * tgxx );
            dt12_yy[ig] = -( v * dy_y12 * tgy - y12 * dv * tg * tgyy );
            dt12_zz[ig] = -( v * dz_y12 * tgz - y12 * dv * tg * tgzz );
1204
            dt12_xy[ig] = -( 0.5 * v * ( dx_y12 * tgy + dy_y12 * tgx )
1205
              - y12 * dv * tg * tgxy );
1206
            dt12_yz[ig] = -( 0.5 * v * ( dy_y12 * tgz + dz_y12 * tgy )
1207
              - y12 * dv * tg * tgyz );
1208
            dt12_xz[ig] = -( 0.5 * v * ( dx_y12 * tgz + dz_y12 * tgx )
1209 1210 1211 1212 1213
              - y12 * dv * tg * tgxz );

            dt13_xx[ig] = -( v * dx_y13 * tgx - y13 * dv * tg * tgxx );
            dt13_yy[ig] = -( v * dy_y13 * tgy - y13 * dv * tg * tgyy );
            dt13_zz[ig] = -( v * dz_y13 * tgz - y13 * dv * tg * tgzz );
1214
            dt13_xy[ig] = -( 0.5 * v * ( dx_y13 * tgy + dy_y13 * tgx )
1215
              - y13 * dv * tg * tgxy );
1216
            dt13_yz[ig] = -( 0.5 * v * ( dy_y13 * tgz + dz_y13 * tgy )
1217
              - y13 * dv * tg * tgyz );
1218
            dt13_xz[ig] = -( 0.5 * v * ( dx_y13 * tgz + dz_y13 * tgx )
1219 1220 1221 1222 1223
              - y13 * dv * tg * tgxz );

            dt14_xx[ig] = -( v * dx_y14 * tgx - y14 * dv * tg * tgxx );
            dt14_yy[ig] = -( v * dy_y14 * tgy - y14 * dv * tg * tgyy );
            dt14_zz[ig] = -( v * dz_y14 * tgz - y14 * dv * tg * tgzz );
1224
            dt14_xy[ig] = -( 0.5 * v * ( dx_y14 * tgy + dy_y14 * tgx )
1225
              - y14 * dv * tg * tgxy );
1226
            dt14_yz[ig] = -( 0.5 * v * ( dy_y14 * tgz + dz_y14 * tgy )
1227
              - y14 * dv * tg * tgyz );
1228
            dt14_xz[ig] = -( 0.5 * v * ( dx_y14 * tgz + dz_y14 * tgx )
1229 1230 1231 1232 1233
              - y14 * dv * tg * tgxz );

            dt15_xx[ig] = -( v * dx_y15 * tgx - y15 * dv * tg * tgxx );
            dt15_yy[ig] = -( v * dy_y15 * tgy - y15 * dv * tg * tgyy );
            dt15_zz[ig] = -( v * dz_y15 * tgz - y15 * dv * tg * tgzz );
1234
            dt15_xy[ig] = -( 0.5 * v * ( dx_y15 * tgy + dy_y15 * tgx )
1235
              - y15 * dv * tg * tgxy );
1236
            dt15_yz[ig] = -( 0.5 * v * ( dy_y15 * tgz + dz_y15 * tgy )
1237
              - y15 * dv * tg * tgyz );
1238
            dt15_xz[ig] = -( 0.5 * v * ( dx_y15 * tgz + dz_y15 * tgx )
1239 1240 1241 1242 1243 1244 1245
              - y15 * dv * tg * tgxz );

          }

          ilm += 2 * l + 1;

        } // l == 3
Francois Gygi committed
1246 1247 1248 1249 1250 1251 1252 1253
        else
        {
          assert(false);
        }
      } // l != lloc[is]
    } // for l
    assert(ilm == s->nlm());
  }
1254
  tmap["update_twnl"].stop();
Francois Gygi committed
1255 1256 1257
}

////////////////////////////////////////////////////////////////////////////////
1258
double NonLocalPotential::energy(bool compute_hpsi, SlaterDet& dsd,
1259
    bool compute_forces, vector<vector<double> >& fion_enl,
1260
    bool compute_stress, valarray<double>& sigma_enl)
Francois Gygi committed
1261
{
1262 1263 1264 1265 1266 1267 1268 1269 1270
  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;
  sigma_enl = 0.0;

  if ( nspnl == 0 ) return 0.0;

  double enl = 0.0;
  double tsum[6] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 };
1271
  const vector<double>& occ = sd_.occ();
1272 1273 1274 1275 1276 1277 1278 1279
  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
1280 1281 1282 1283

  const double omega = basis_.cell().volume();
  assert(omega != 0.0);
  const double omega_inv = 1.0 / omega;
1284 1285 1286

  for ( int is = 0; is < nsp; is++ )
  {
1287
    Species *s = atoms_.species_list[is];
Francois Gygi committed
1288 1289
    if ( npr[is] > 0 ) // species is is non-local
    {
1290 1291 1292 1293 1294 1295 1296 1297 1298
      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
1299 1300 1301 1302 1303
      // 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);
1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320
      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
1321
        double* kpgx = const_cast<double*>(basis_.kpgx_ptr(0));
1322
        dgemm(&cn,&cn,(int*)&ngwl,(int*)&ia_block_size,&k,&mone,
Francois Gygi committed
1323
              kpgx,(int*)&ngwl, &tau[is][3*iastart],&k,
1324 1325 1326
              &zero,&gr[0],(int*)&ngwl);

        int len = ia_block_size * ngwl;
1327
#if USE_MASSV
1328
        vsincos(&sgr[0],&cgr[0],&gr[0],&len);
1329
#else
1330 1331 1332 1333 1334 1335
        for ( int i = 0; i < len; i++ )
        {
          const double arg = gr[i];
          sgr[i] = sin(arg);
          cgr[i] = cos(arg);
        }
1336
#endif
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 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380
        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