ReferenceLJCoulombIxn.cpp 29.7 KB
Newer Older
1

2
/* Portions copyright (c) 2006-2013 Stanford University and Simbios.
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
 * Contributors: Pande Group
 *
 * Permission is hereby granted, free of charge, to any person obtaining
 * a copy of this software and associated documentation files (the
 * "Software"), to deal in the Software without restriction, including
 * without limitation the rights to use, copy, modify, merge, publish,
 * distribute, sublicense, and/or sell copies of the Software, and to
 * permit persons to whom the Software is furnished to do so, subject
 * to the following conditions:
 *
 * The above copyright notice and this permission notice shall be included
 * in all copies or substantial portions of the Software.
 *
 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
 * OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
 * MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
 * IN NO EVENT SHALL THE AUTHORS, CONTRIBUTORS OR COPYRIGHT HOLDERS BE
 * LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION
 * OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION
 * WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
 */

#include <string.h>
#include <sstream>
27
#include <complex>
28
#include <algorithm>
29
#include <iostream>
30

31
#include "SimTKOpenMMUtilities.h"
32
33
#include "ReferenceLJCoulombIxn.h"
#include "ReferenceForce.h"
34
#include "ReferencePME.h"
peastman's avatar
peastman committed
35
#include "openmm/OpenMMException.h"
36

37
38
// In case we're using some primitive version of Visual Studio this will
// make sure that erf() and erfc() are defined.
39
#include "openmm/internal/MSVC_erfc.h"
40

41
using std::set;
42
using std::vector;
43
using namespace OpenMM;
44

45
46
47
48
49
50
/**---------------------------------------------------------------------------------------

   ReferenceLJCoulombIxn constructor

   --------------------------------------------------------------------------------------- */

51
ReferenceLJCoulombIxn::ReferenceLJCoulombIxn() : cutoff(false), useSwitch(false), periodic(false), ewald(false), pme(false), ljpme(false) {
52

53
    // ---------------------------------------------------------------------------------------
54

55
    // static const char* methodName = "\nReferenceLJCoulombIxn::ReferenceLJCoulombIxn";
56

57
    // ---------------------------------------------------------------------------------------
58
59
60
61
62
63
64
65
66

}

/**---------------------------------------------------------------------------------------

   ReferenceLJCoulombIxn destructor

   --------------------------------------------------------------------------------------- */

67
ReferenceLJCoulombIxn::~ReferenceLJCoulombIxn() {
68

69
    // ---------------------------------------------------------------------------------------
70

71
    // static const char* methodName = "\nReferenceLJCoulombIxn::~ReferenceLJCoulombIxn";
72

73
    // ---------------------------------------------------------------------------------------
74
75
76

}

77
/**---------------------------------------------------------------------------------------
78
79
80
81
82
83
84
85
86

     Set the force to use a cutoff.

     @param distance            the cutoff distance
     @param neighbors           the neighbor list to use
     @param solventDielectric   the dielectric constant of the bulk solvent

     --------------------------------------------------------------------------------------- */

87
void ReferenceLJCoulombIxn::setUseCutoff(RealOpenMM distance, const OpenMM::NeighborList& neighbors, RealOpenMM solventDielectric) {
88

89
90
91
    cutoff = true;
    cutoffDistance = distance;
    neighborList = &neighbors;
92
93
    krf = pow(cutoffDistance, -3.0)*(solventDielectric-1.0)/(2.0*solventDielectric+1.0);
    crf = (1.0/cutoffDistance)*(3.0*solventDielectric)/(2.0*solventDielectric+1.0);
94
}
95

96
97
98
99
100
101
102
103
/**---------------------------------------------------------------------------------------

   Set the force to use a switching function on the Lennard-Jones interaction.

   @param distance            the switching distance

   --------------------------------------------------------------------------------------- */

104
void ReferenceLJCoulombIxn::setUseSwitchingFunction(RealOpenMM distance) {
105
106
107
108
    useSwitch = true;
    switchingDistance = distance;
}

109
/**---------------------------------------------------------------------------------------
110
111
112
113
114

     Set the force to use periodic boundary conditions.  This requires that a cutoff has
     also been set, and the smallest side of the periodic box is at least twice the cutoff
     distance.

115
     @param vectors    the vectors defining the periodic box
116
117
118

     --------------------------------------------------------------------------------------- */

119
void ReferenceLJCoulombIxn::setPeriodic(OpenMM::RealVec* vectors) {
120
121

    assert(cutoff);
122
123
124
    assert(vectors[0][0] >= 2.0*cutoffDistance);
    assert(vectors[1][1] >= 2.0*cutoffDistance);
    assert(vectors[2][2] >= 2.0*cutoffDistance);
125
    periodic = true;
126
127
128
    periodicBoxVectors[0] = vectors[0];
    periodicBoxVectors[1] = vectors[1];
    periodicBoxVectors[2] = vectors[2];
129
}
130

131
/**---------------------------------------------------------------------------------------
132
133
134
135
136
137
138
139
140
141

     Set the force to use Ewald summation.

     @param alpha  the Ewald separation parameter
     @param kmaxx  the largest wave vector in the x direction
     @param kmaxy  the largest wave vector in the y direction
     @param kmaxz  the largest wave vector in the z direction

     --------------------------------------------------------------------------------------- */

142
143
144
145
146
147
148
void ReferenceLJCoulombIxn::setUseEwald(RealOpenMM alpha, int kmaxx, int kmaxy, int kmaxz) {
    alphaEwald = alpha;
    numRx = kmaxx;
    numRy = kmaxy;
    numRz = kmaxz;
    ewald = true;
}
149

150
/**---------------------------------------------------------------------------------------
151
152
153
154

     Set the force to use Particle-Mesh Ewald (PME) summation.

     @param alpha  the Ewald separation parameter
155
     @param gridSize the dimensions of the mesh
156
157
158

     --------------------------------------------------------------------------------------- */

159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
void ReferenceLJCoulombIxn::setUsePME(RealOpenMM alpha, int meshSize[3]) {
    alphaEwald = alpha;
    meshDim[0] = meshSize[0];
    meshDim[1] = meshSize[1];
    meshDim[2] = meshSize[2];
    pme = true;
}

/**---------------------------------------------------------------------------------------

     Set the force to use Particle-Mesh Ewald (PME) summation for dispersion and terms.

     @param dalpha  the dispersion Ewald separation parameter
     @param dgridSize the dimensions of the dispersion mesh

     --------------------------------------------------------------------------------------- */

void ReferenceLJCoulombIxn::setUseLJPME(RealOpenMM dalpha, int dmeshSize[3]) {
    alphaDispersionEwald = dalpha;
    dispersionMeshDim[0] = dmeshSize[0];
    dispersionMeshDim[1] = dmeshSize[1];
    dispersionMeshDim[2] = dmeshSize[2];
    ljpme = true;
}
183

184
185
186
187
188
189
190
/**---------------------------------------------------------------------------------------

   Calculate Ewald ixn

   @param numberOfAtoms    number of atoms
   @param atomCoordinates  atom coordinates
   @param atomParameters   atom parameters                             atomParameters[atomIndex][paramterIndex]
191
192
   @param exclusions       atom exclusion indices
                           exclusions[atomIndex] contains the list of exclusions for that atom
193
194
195
196
   @param fixedParameters  non atom parameters (not currently used)
   @param forces           force array (forces added)
   @param energyByAtom     atom energy
   @param totalEnergy      total energy
197
198
   @param includeDirect      true if direct space interactions should be included
   @param includeReciprocal  true if reciprocal space interactions should be included
199
200

   --------------------------------------------------------------------------------------- */
201

202
void ReferenceLJCoulombIxn::calculateEwaldIxn(int numberOfAtoms, vector<RealVec>& atomCoordinates,
203
204
205
                                              RealOpenMM** atomParameters, vector<set<int> >& exclusions,
                                              RealOpenMM* fixedParameters, vector<RealVec>& forces,
                                              RealOpenMM* energyByAtom, RealOpenMM* totalEnergy, bool includeDirect, bool includeReciprocal) const {
206
207
    typedef std::complex<RealOpenMM> d_complex;

208
    static const RealOpenMM epsilon     =  1.0;
209
    static const RealOpenMM one         =  1.0;
210
211
    static const RealOpenMM six         =  6.0;
    static const RealOpenMM twelve      = 12.0;
212

213
    int kmax                            = (ewald ? std::max(numRx, std::max(numRy,numRz)) : 0);
214
    RealOpenMM  factorEwald             = -1 / (4*alphaEwald*alphaEwald);
Peter Eastman's avatar
Peter Eastman committed
215
216
    RealOpenMM SQRT_PI                  = sqrt(PI_M);
    RealOpenMM TWO_PI                   = 2.0 * PI_M;
217
    RealOpenMM recipCoeff               = (RealOpenMM)(ONE_4PI_EPS0*4*PI_M/(periodicBoxVectors[0][0] * periodicBoxVectors[1][1] * periodicBoxVectors[2][2]) /epsilon);
218

219
    RealOpenMM totalSelfEwaldEnergy     = 0.0;
220
221
    RealOpenMM realSpaceEwaldEnergy     = 0.0;
    RealOpenMM recipEnergy              = 0.0;
222
    RealOpenMM recipDispersionEnergy    = 0.0;
223
    RealOpenMM totalRecipEnergy         = 0.0;
224
    RealOpenMM vdwEnergy                = 0.0;
225

226
227
228
229
230
231
232
233
234
    // A couple of sanity checks for
    if(ljpme && useSwitch)
        throw OpenMMException("Switching cannot be used with LJPME");
    if(ljpme && !pme)
        throw OpenMMException("LJPME has been set, without PME being set");

    // **************************************************************************************
    // SELF ENERGY
    // **************************************************************************************
235

236
    if (includeReciprocal) {
237
        for (int atomID = 0; atomID < numberOfAtoms; atomID++) {
238
            RealOpenMM selfEwaldEnergy       = (RealOpenMM) (ONE_4PI_EPS0*atomParameters[atomID][QIndex]*atomParameters[atomID][QIndex] * alphaEwald/SQRT_PI);
239
240
241
242
            if(ljpme) {
                // Dispersion self term
                selfEwaldEnergy -= pow(alphaDispersionEwald, 6.0) * 64.0*pow(atomParameters[atomID][SigIndex], 6.0) * pow(atomParameters[atomID][EpsIndex], 2.0) / 12.0;
            }
243
            totalSelfEwaldEnergy            -= selfEwaldEnergy;
244
            if (energyByAtom) {
245
246
                energyByAtom[atomID]        -= selfEwaldEnergy;
            }
247
        }
248
249
    }

250
    if (totalEnergy) {
251
252
253
        *totalEnergy += totalSelfEwaldEnergy;
    }

254
255
256
    // **************************************************************************************
    // RECIPROCAL SPACE EWALD ENERGY AND FORCES
    // **************************************************************************************
257
    // PME
258

259
260
    if (pme && includeReciprocal) {
        pme_t          pmedata; /* abstract handle for PME data */
261

262
        pme_init(&pmedata,alphaEwald,numberOfAtoms,meshDim,5,1);
263

264
265
266
267
        vector<RealOpenMM> charges(numberOfAtoms);
        for (int i = 0; i < numberOfAtoms; i++)
            charges[i] = atomParameters[i][QIndex];
        pme_exec(pmedata,atomCoordinates,forces,charges,periodicBoxVectors,&recipEnergy);
268

269
270
        if (totalEnergy)
            *totalEnergy += recipEnergy;
271

272
273
274
        if (energyByAtom)
            for (int n = 0; n < numberOfAtoms; n++)
                energyByAtom[n] += recipEnergy;
275

276
        pme_destroy(pmedata);
277

278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
        if (ljpme) {
            // Dispersion reciprocal space terms
            pme_init(&pmedata,alphaDispersionEwald,numberOfAtoms,dispersionMeshDim,5,1);

            std::vector<RealVec> dpmeforces;
            for (int i = 0; i < numberOfAtoms; i++){
                charges[i] = 8.0*pow(atomParameters[i][SigIndex], 3.0) * atomParameters[i][EpsIndex];
                dpmeforces.push_back(RealVec());
            }
            pme_exec_dpme(pmedata,atomCoordinates,dpmeforces,charges,periodicBoxVectors,&recipDispersionEnergy);
            for (int i = 0; i < numberOfAtoms; i++){
                forces[i][0] -= 2.0*dpmeforces[i][0];
                forces[i][1] -= 2.0*dpmeforces[i][1];
                forces[i][2] -= 2.0*dpmeforces[i][2];
            }
            if (totalEnergy)
                *totalEnergy += recipDispersionEnergy;

            if (energyByAtom)
                for (int n = 0; n < numberOfAtoms; n++)
                    energyByAtom[n] += recipDispersionEnergy;
            pme_destroy(pmedata);
        }
    }
302
303
    // Ewald method

304
    else if (ewald && includeReciprocal) {
305

306
        // setup reciprocal box
307

308
        RealOpenMM recipBoxSize[3] = { TWO_PI / periodicBoxVectors[0][0], TWO_PI / periodicBoxVectors[1][1], TWO_PI / periodicBoxVectors[2][2]};
309
310


311
        // setup K-vectors
312

313
314
315
316
#define EIR(x, y, z) eir[(x)*numberOfAtoms*3+(y)*3+z]
        vector<d_complex> eir(kmax*numberOfAtoms*3);
        vector<d_complex> tab_xy(numberOfAtoms);
        vector<d_complex> tab_qxyz(numberOfAtoms);
317

318
319
        if (kmax < 1)
            throw OpenMMException("kmax for Ewald summation < 1");
320

321
322
323
        for (int i = 0; (i < numberOfAtoms); i++) {
            for (int m = 0; (m < 3); m++)
                EIR(0, i, m) = d_complex(1,0);
324

325
326
327
            for (int m=0; (m<3); m++)
                EIR(1, i, m) = d_complex(cos(atomCoordinates[i][m]*recipBoxSize[m]),
                                         sin(atomCoordinates[i][m]*recipBoxSize[m]));
328

329
330
331
332
            for (int j=2; (j<kmax); j++)
                for (int m=0; (m<3); m++)
                    EIR(j, i, m) = EIR(j-1, i, m) * EIR(1, i, m);
        }
333

334
        // calculate reciprocal space energy and forces
335

336
337
        int lowry = 0;
        int lowrz = 1;
338

339
        for (int rx = 0; rx < numRx; rx++) {
340

341
            RealOpenMM kx = rx * recipBoxSize[0];
342

343
            for (int ry = lowry; ry < numRy; ry++) {
344

345
                RealOpenMM ky = ry * recipBoxSize[1];
346

347
348
349
350
                if (ry >= 0) {
                    for (int n = 0; n < numberOfAtoms; n++)
                        tab_xy[n] = EIR(rx, n, 0) * EIR(ry, n, 1);
                }
351

352
353
354
355
                else {
                    for (int n = 0; n < numberOfAtoms; n++)
                        tab_xy[n]= EIR(rx, n, 0) * conj (EIR(-ry, n, 1));
                }
356

357
                for (int rz = lowrz; rz < numRz; rz++) {
358

359
360
361
362
                    if (rz >= 0) {
                        for (int n = 0; n < numberOfAtoms; n++)
                            tab_qxyz[n] = atomParameters[n][QIndex] * (tab_xy[n] * EIR(rz, n, 2));
                    }
363

364
365
366
367
                    else {
                        for (int n = 0; n < numberOfAtoms; n++)
                            tab_qxyz[n] = atomParameters[n][QIndex] * (tab_xy[n] * conj(EIR(-rz, n, 2)));
                    }
368

369
370
                    RealOpenMM cs = 0.0f;
                    RealOpenMM ss = 0.0f;
371

372
373
374
375
                    for (int n = 0; n < numberOfAtoms; n++) {
                        cs += tab_qxyz[n].real();
                        ss += tab_qxyz[n].imag();
                    }
376

377
378
379
                    RealOpenMM kz = rz * recipBoxSize[2];
                    RealOpenMM k2 = kx * kx + ky * ky + kz * kz;
                    RealOpenMM ak = exp(k2*factorEwald) / k2;
380

381
382
383
384
385
386
                    for (int n = 0; n < numberOfAtoms; n++) {
                        RealOpenMM force = ak * (cs * tab_qxyz[n].imag() - ss * tab_qxyz[n].real());
                        forces[n][0] += 2 * recipCoeff * force * kx ;
                        forces[n][1] += 2 * recipCoeff * force * ky ;
                        forces[n][2] += 2 * recipCoeff * force * kz ;
                    }
387

388
389
                    recipEnergy       = recipCoeff * ak * (cs * cs + ss * ss);
                    totalRecipEnergy += recipEnergy;
390

391
392
                    if (totalEnergy)
                        *totalEnergy += recipEnergy;
393

394
395
396
                    if (energyByAtom)
                        for (int n = 0; n < numberOfAtoms; n++)
                            energyByAtom[n] += recipEnergy;
397

398
399
400
401
                    lowrz = 1 - numRz;
                }
                lowry = 1 - numRy;
            }
402
403
404
        }
    }

405
406
407
    // **************************************************************************************
    // SHORT-RANGE ENERGY AND FORCES
    // **************************************************************************************
408

409
410
    if (!includeDirect)
        return;
411
412
413
    RealOpenMM totalVdwEnergy            = 0.0f;
    RealOpenMM totalRealSpaceEwaldEnergy = 0.0f;

414

415
    for (int i = 0; i < (int) neighborList->size(); i++) {
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
        OpenMM::AtomPair pair = (*neighborList)[i];
        int ii = pair.first;
        int jj = pair.second;

        RealOpenMM deltaR[2][ReferenceForce::LastDeltaRIndex];
        ReferenceForce::getDeltaRPeriodic(atomCoordinates[jj], atomCoordinates[ii], periodicBoxVectors, deltaR[0]);
        RealOpenMM r         = deltaR[0][ReferenceForce::RIndex];
        RealOpenMM inverseR  = one/(deltaR[0][ReferenceForce::RIndex]);
        RealOpenMM switchValue = 1, switchDeriv = 0;
        if (useSwitch && r > switchingDistance) {
            RealOpenMM t = (r-switchingDistance)/(cutoffDistance-switchingDistance);
            switchValue = 1+t*t*t*(-10+t*(15-t*6));
            switchDeriv = t*t*(-30+t*(60-t*30))/(cutoffDistance-switchingDistance);
        }
        RealOpenMM alphaR    = alphaEwald * r;


        RealOpenMM dEdR      = (RealOpenMM) (ONE_4PI_EPS0 * atomParameters[ii][QIndex] * atomParameters[jj][QIndex] * inverseR * inverseR * inverseR);
        dEdR      = (RealOpenMM) (dEdR * (erfc(alphaR) + 2 * alphaR * exp (- alphaR * alphaR) / SQRT_PI));

        RealOpenMM sig       = atomParameters[ii][SigIndex] +  atomParameters[jj][SigIndex];
        RealOpenMM sig2      = inverseR*sig;
        sig2     *= sig2;
        RealOpenMM sig6      = sig2*sig2*sig2;
        RealOpenMM eps       = atomParameters[ii][EpsIndex]*atomParameters[jj][EpsIndex];
        dEdR     += switchValue*eps*(twelve*sig6 - six)*sig6*inverseR*inverseR;
        vdwEnergy = eps*(sig6-one)*sig6;

        if (ljpme) {
            RealOpenMM dalphaR   = alphaDispersionEwald * r;
            RealOpenMM dar2 = dalphaR*dalphaR;
            RealOpenMM dar4 = dar2*dar2;
            RealOpenMM dar6 = dar4*dar2;
            RealOpenMM inverseR2 = inverseR*inverseR;
            RealOpenMM c6i = 8.0*pow(atomParameters[ii][SigIndex], 3.0) * atomParameters[ii][EpsIndex];
            RealOpenMM c6j = 8.0*pow(atomParameters[jj][SigIndex], 3.0) * atomParameters[jj][EpsIndex];
            // For the energies and forces, we first add the regular Lorentz−Berthelot terms.  The C12 term is treated as usual
            // but we then subtract out (remembering that the C6 term is negative) the multiplicative C6 term that has been
            // computed in real space.  Finally, we add a potential shift term to account for the difference between the LB
            // and multiplicative functional forms at the cutoff.
            RealOpenMM emult = c6i*c6j*inverseR2*inverseR2*inverseR2*(1.0 - EXP(-dar2) * (1.0 + dar2 + 0.5*dar4));
            dEdR += 6.0*c6i*c6j*inverseR2*inverseR2*inverseR2*inverseR2*(1.0 - EXP(-dar2) * (1.0 + dar2 + 0.5*dar4 + dar6/6.0));
            // The additive part of the potential shift
            RealOpenMM inverseCut2 = 1.0/(cutoffDistance*cutoffDistance);
            RealOpenMM inverseCut6 = inverseCut2*inverseCut2*inverseCut2;
            sig2 = atomParameters[ii][SigIndex] +  atomParameters[jj][SigIndex];
            sig2 *= sig2;
            sig6 = sig2*sig2*sig2;
            RealOpenMM potentialshift = eps*sig6*inverseCut6;
            dalphaR   = alphaDispersionEwald * cutoffDistance;
            dar2 = dalphaR*dalphaR;
            dar4 = dar2*dar2;
            // The multiplicative part of the potential shift
            potentialshift -= c6i*c6j*inverseCut6*(1.0 - EXP(-dar2) * (1.0 + dar2 + 0.5*dar4));
            vdwEnergy += emult + potentialshift;
        }

        if (useSwitch) {
            dEdR -= vdwEnergy*switchDeriv*inverseR;
            vdwEnergy *= switchValue;
        }

        // accumulate forces

        for (int kk = 0; kk < 3; kk++) {
            RealOpenMM force  = dEdR*deltaR[0][kk];
            forces[ii][kk]   += force;
            forces[jj][kk]   -= force;
        }

        // accumulate energies

        realSpaceEwaldEnergy        = (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*inverseR*erfc(alphaR));

        totalVdwEnergy             += vdwEnergy;
        totalRealSpaceEwaldEnergy  += realSpaceEwaldEnergy;
492

493
        if (energyByAtom) {
494
495
            energyByAtom[ii] += realSpaceEwaldEnergy + vdwEnergy;
            energyByAtom[jj] += realSpaceEwaldEnergy + vdwEnergy;
496
        }
497

498
    }
499

500
    if (totalEnergy)
501
502
        *totalEnergy += totalRealSpaceEwaldEnergy + totalVdwEnergy;

Peter Eastman's avatar
Peter Eastman committed
503
504
    // Now subtract off the exclusions, since they were implicitly included in the reciprocal space sum.

505
    RealOpenMM totalExclusionEnergy = 0.0f;
506
    const double TWO_OVER_SQRT_PI = 2/sqrt(PI_M);
Peter Eastman's avatar
Peter Eastman committed
507
    for (int i = 0; i < numberOfAtoms; i++)
508
509
        for (set<int>::const_iterator iter = exclusions[i].begin(); iter != exclusions[i].end(); ++iter) {
            if (*iter > i) {
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
                int ii = i;
                int jj = *iter;

                RealOpenMM deltaR[2][ReferenceForce::LastDeltaRIndex];
                ReferenceForce::getDeltaR(atomCoordinates[jj], atomCoordinates[ii], deltaR[0]);
                RealOpenMM r         = deltaR[0][ReferenceForce::RIndex];
                RealOpenMM inverseR  = one/(deltaR[0][ReferenceForce::RIndex]);
                RealOpenMM alphaR    = alphaEwald * r;
                if (erf(alphaR) > 1e-6) {
                    RealOpenMM dEdR      = (RealOpenMM) (ONE_4PI_EPS0 * atomParameters[ii][QIndex] * atomParameters[jj][QIndex] * inverseR * inverseR * inverseR);
                    dEdR      = (RealOpenMM) (dEdR * (erf(alphaR) - 2 * alphaR * exp (- alphaR * alphaR) / SQRT_PI));

                    // accumulate forces

                    for (int kk = 0; kk < 3; kk++) {
                        RealOpenMM force  = dEdR*deltaR[0][kk];
                        forces[ii][kk]   -= force;
                        forces[jj][kk]   += force;
                    }

                    // accumulate energies

                    realSpaceEwaldEnergy = (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*inverseR*erf(alphaR));
                }
                else {
                    realSpaceEwaldEnergy = (RealOpenMM) (alphaEwald*TWO_OVER_SQRT_PI*ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]);
                }

                if(ljpme){
                    // Dispersion terms.  Here we just back out the reciprocal space terms, and don't add any extra real space terms.
                    RealOpenMM dalphaR   = alphaDispersionEwald * r;
                    RealOpenMM inverseR2 = inverseR*inverseR;
                    RealOpenMM dar2 = dalphaR*dalphaR;
                    RealOpenMM dar4 = dar2*dar2;
                    RealOpenMM dar6 = dar4*dar2;
                    RealOpenMM c6i = 8.0*pow(atomParameters[ii][SigIndex], 3.0) * atomParameters[ii][EpsIndex];
                    RealOpenMM c6j = 8.0*pow(atomParameters[jj][SigIndex], 3.0) * atomParameters[jj][EpsIndex];
                    realSpaceEwaldEnergy -= c6i*c6j*inverseR2*inverseR2*inverseR2*(1.0 - EXP(-dar2) * (1.0 + dar2 + 0.5*dar4));
                    RealOpenMM dEdR = -6.0*c6i*c6j*inverseR2*inverseR2*inverseR2*inverseR2*(1.0 - EXP(-dar2) * (1.0 + dar2 + 0.5*dar4 + dar6/6.0));
                    for (int kk = 0; kk < 3; kk++) {
                        RealOpenMM force  = dEdR*deltaR[0][kk];
                        forces[ii][kk]   -= force;
                        forces[jj][kk]   += force;
                    }
                }

                totalExclusionEnergy += realSpaceEwaldEnergy;
                if (energyByAtom) {
                    energyByAtom[ii] -= realSpaceEwaldEnergy;
                    energyByAtom[jj] -= realSpaceEwaldEnergy;
                }
Peter Eastman's avatar
Peter Eastman committed
561
            }
562
        }
Peter Eastman's avatar
Peter Eastman committed
563

564
    if (totalEnergy)
565
        *totalEnergy -= totalExclusionEnergy;
566
567
568
}


569
570
571
572
573
574
575
/**---------------------------------------------------------------------------------------

   Calculate LJ Coulomb pair ixn

   @param numberOfAtoms    number of atoms
   @param atomCoordinates  atom coordinates
   @param atomParameters   atom parameters                             atomParameters[atomIndex][paramterIndex]
576
577
   @param exclusions       atom exclusion indices
                           exclusions[atomIndex] contains the list of exclusions for that atom
578
579
580
581
   @param fixedParameters  non atom parameters (not currently used)
   @param forces           force array (forces added)
   @param energyByAtom     atom energy
   @param totalEnergy      total energy
582
583
   @param includeDirect      true if direct space interactions should be included
   @param includeReciprocal  true if reciprocal space interactions should be included
584
585

   --------------------------------------------------------------------------------------- */
586

587
void ReferenceLJCoulombIxn::calculatePairIxn(int numberOfAtoms, vector<RealVec>& atomCoordinates,
588
                                             RealOpenMM** atomParameters, vector<set<int> >& exclusions,
589
                                             RealOpenMM* fixedParameters, vector<RealVec>& forces,
590
                                             RealOpenMM* energyByAtom, RealOpenMM* totalEnergy, bool includeDirect, bool includeReciprocal) const {
591

592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
    if (ewald || pme || ljpme) {
        calculateEwaldIxn(numberOfAtoms, atomCoordinates, atomParameters, exclusions, fixedParameters, forces, energyByAtom,
                          totalEnergy, includeDirect, includeReciprocal);
        return;
    }
    if (!includeDirect)
        return;
    if (cutoff) {
        for (int i = 0; i < (int) neighborList->size(); i++) {
            OpenMM::AtomPair pair = (*neighborList)[i];
            calculateOneIxn(pair.first, pair.second, atomCoordinates, atomParameters, forces, energyByAtom, totalEnergy);
        }
    }
    else {
        for (int ii = 0; ii < numberOfAtoms; ii++) {
            // loop over atom pairs

            for (int jj = ii+1; jj < numberOfAtoms; jj++)
                if (exclusions[jj].find(ii) == exclusions[jj].end())
                    calculateOneIxn(ii, jj, atomCoordinates, atomParameters, forces, energyByAtom, totalEnergy);
        }
    }
614
}
615

616
/**---------------------------------------------------------------------------------------
617

618
     Calculate LJ Coulomb pair ixn between two atoms
619

620
621
622
623
624
625
626
     @param ii               the index of the first atom
     @param jj               the index of the second atom
     @param atomCoordinates  atom coordinates
     @param atomParameters   atom parameters (charges, c6, c12, ...)     atomParameters[atomIndex][paramterIndex]
     @param forces           force array (forces added)
     @param energyByAtom     atom energy
     @param totalEnergy      total energy
627

628
     --------------------------------------------------------------------------------------- */
629

630
void ReferenceLJCoulombIxn::calculateOneIxn(int ii, int jj, vector<RealVec>& atomCoordinates,
631
632
                                            RealOpenMM** atomParameters, vector<RealVec>& forces,
                                            RealOpenMM* energyByAtom, RealOpenMM* totalEnergy) const {
633

634
    // ---------------------------------------------------------------------------------------
635

636
    static const std::string methodName = "\nReferenceLJCoulombIxn::calculateOneIxn";
637

638
    // ---------------------------------------------------------------------------------------
639

640
    // constants -- reduce Visual Studio warnings regarding conversions between float & double
641

642
643
644
645
646
647
648
    static const RealOpenMM zero        =  0.0;
    static const RealOpenMM one         =  1.0;
    static const RealOpenMM two         =  2.0;
    static const RealOpenMM three       =  3.0;
    static const RealOpenMM six         =  6.0;
    static const RealOpenMM twelve      = 12.0;
    static const RealOpenMM oneM        = -1.0;
649

650
    static const int threeI             = 3;
651

652
653
654
655
656
657
658
    static const int LastAtomIndex      = 2;

    RealOpenMM deltaR[2][ReferenceForce::LastDeltaRIndex];

    // get deltaR, R2, and R between 2 atoms

    if (periodic)
659
        ReferenceForce::getDeltaRPeriodic(atomCoordinates[jj], atomCoordinates[ii], periodicBoxVectors, deltaR[0]);
660
    else
661
        ReferenceForce::getDeltaR(atomCoordinates[jj], atomCoordinates[ii], deltaR[0]);
662
663
664

    RealOpenMM r2        = deltaR[0][ReferenceForce::R2Index];
    RealOpenMM inverseR  = one/(deltaR[0][ReferenceForce::RIndex]);
665
666
667
668
669
670
671
672
673
    RealOpenMM switchValue = 1, switchDeriv = 0;
    if (useSwitch) {
        RealOpenMM r = deltaR[0][ReferenceForce::RIndex];
        if (r > switchingDistance) {
            RealOpenMM t = (r-switchingDistance)/(cutoffDistance-switchingDistance);
            switchValue = 1+t*t*t*(-10+t*(15-t*6));
            switchDeriv = t*t*(-30+t*(60-t*30))/(cutoffDistance-switchingDistance);
        }
    }
674
675
    RealOpenMM sig       = atomParameters[ii][SigIndex] +  atomParameters[jj][SigIndex];
    RealOpenMM sig2      = inverseR*sig;
676
    sig2     *= sig2;
677
678
679
    RealOpenMM sig6      = sig2*sig2*sig2;

    RealOpenMM eps       = atomParameters[ii][EpsIndex]*atomParameters[jj][EpsIndex];
680
    RealOpenMM dEdR      = switchValue*eps*(twelve*sig6 - six)*sig6;
681
682
683
684
685
686
687
688
689
690
691
692
693
694
    if (cutoff)
        dEdR += (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*(inverseR-2.0f*krf*r2));
    else
        dEdR += (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*inverseR);
    dEdR     *= inverseR*inverseR;
    RealOpenMM energy = eps*(sig6-one)*sig6;
    if (useSwitch) {
        dEdR -= energy*switchDeriv*inverseR;
        energy *= switchValue;
    }
    if (cutoff)
        energy += (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*(inverseR+krf*r2-crf));
    else
        energy += (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*inverseR);
695

696
697
    // accumulate forces

698
    for (int kk = 0; kk < 3; kk++) {
699
700
701
        RealOpenMM force  = dEdR*deltaR[0][kk];
        forces[ii][kk]   += force;
        forces[jj][kk]   -= force;
702
703
704
705
    }

    // accumulate energies

706
    if (totalEnergy)
707
        *totalEnergy += energy;
708
    if (energyByAtom) {
709
710
        energyByAtom[ii] += energy;
        energyByAtom[jj] += energy;
711
    }
712
}
713