ReferenceLJCoulombIxn.cpp 24.1 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

29
30
31
#include "../SimTKUtilities/SimTKOpenMMCommon.h"
#include "../SimTKUtilities/SimTKOpenMMLog.h"
#include "../SimTKUtilities/SimTKOpenMMUtilities.h"
32
33
#include "ReferenceLJCoulombIxn.h"
#include "ReferenceForce.h"
34
#include "PME.h"
35

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

40
using std::vector;
41
using OpenMM::RealVec;
42

43
44
45
46
47
48
/**---------------------------------------------------------------------------------------

   ReferenceLJCoulombIxn constructor

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

49
ReferenceLJCoulombIxn::ReferenceLJCoulombIxn( ) : cutoff(false), useSwitch(false), periodic(false), ewald(false), pme(false) {
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74

   // ---------------------------------------------------------------------------------------

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

   // ---------------------------------------------------------------------------------------

}

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

   ReferenceLJCoulombIxn destructor

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

ReferenceLJCoulombIxn::~ReferenceLJCoulombIxn( ){

   // ---------------------------------------------------------------------------------------

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

   // ---------------------------------------------------------------------------------------

}

75
76
77
78
79
80
81
82
83
84
  /**---------------------------------------------------------------------------------------

     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

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

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

87
88
89
    cutoff = true;
    cutoffDistance = distance;
    neighborList = &neighbors;
90
91
    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);
92
93
  }

94
95
96
97
98
99
100
101
102
103
104
105
106
/**---------------------------------------------------------------------------------------

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

   @param distance            the switching distance

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

void ReferenceLJCoulombIxn::setUseSwitchingFunction( RealOpenMM distance ) {
    useSwitch = true;
    switchingDistance = distance;
}

107
108
109
110
111
112
113
114
115
116
  /**---------------------------------------------------------------------------------------

     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.

     @param boxSize             the X, Y, and Z widths of the periodic box

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

117
  void ReferenceLJCoulombIxn::setPeriodic( RealVec& boxSize ) {
118
119
120
121
122
123
124
125
126
127
128

    assert(cutoff);
    assert(boxSize[0] >= 2.0*cutoffDistance);
    assert(boxSize[1] >= 2.0*cutoffDistance);
    assert(boxSize[2] >= 2.0*cutoffDistance);
    periodic = true;
    periodicBoxSize[0] = boxSize[0];
    periodicBoxSize[1] = boxSize[1];
    periodicBoxSize[2] = boxSize[2];
  }

129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
  /**---------------------------------------------------------------------------------------

     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

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

  void ReferenceLJCoulombIxn::setUseEwald(RealOpenMM alpha, int kmaxx, int kmaxy, int kmaxz) {
      alphaEwald = alpha;
      numRx = kmaxx;
      numRy = kmaxy;
      numRz = kmaxz;
      ewald = true;
  }

148
149
150
151
152
  /**---------------------------------------------------------------------------------------

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

     @param alpha  the Ewald separation parameter
153
     @param gridSize the dimensions of the mesh
154
155
156

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

157
  void ReferenceLJCoulombIxn::setUsePME(RealOpenMM alpha, int meshSize[3]) {
158
      alphaEwald = alpha;
159
160
161
      meshDim[0] = meshSize[0];
      meshDim[1] = meshSize[1];
      meshDim[2] = meshSize[2];
162
163
164
      pme = true;
  }

165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
/**---------------------------------------------------------------------------------------

   Calculate Ewald ixn

   @param numberOfAtoms    number of atoms
   @param atomCoordinates  atom coordinates
   @param atomParameters   atom parameters                             atomParameters[atomIndex][paramterIndex]
   @param exclusions       atom exclusion indices                      exclusions[atomIndex][atomToExcludeIndex]
                           exclusions[atomIndex][0] = number of exclusions
                           exclusions[atomIndex][1-no.] = atom indices of atoms to excluded from
                           interacting w/ atom atomIndex
   @param fixedParameters  non atom parameters (not currently used)
   @param forces           force array (forces added)
   @param energyByAtom     atom energy
   @param totalEnergy      total energy
180
181
   @param includeDirect      true if direct space interactions should be included
   @param includeReciprocal  true if reciprocal space interactions should be included
182
183

   --------------------------------------------------------------------------------------- */
184

185
void ReferenceLJCoulombIxn::calculateEwaldIxn(int numberOfAtoms, vector<RealVec>& atomCoordinates,
186
                                             RealOpenMM** atomParameters, int** exclusions,
187
                                             RealOpenMM* fixedParameters, vector<RealVec>& forces,
188
                                             RealOpenMM* energyByAtom, RealOpenMM* totalEnergy, bool includeDirect, bool includeReciprocal) const {
189
190
    typedef std::complex<RealOpenMM> d_complex;

191
    static const RealOpenMM epsilon     =  1.0;
192
    static const RealOpenMM one         =  1.0;
193
194
    static const RealOpenMM six         =  6.0;
    static const RealOpenMM twelve      = 12.0;
195

196
    int kmax                            = (ewald ? std::max(numRx, std::max(numRy,numRz)) : 0);
197
    RealOpenMM  factorEwald             = -1 / (4*alphaEwald*alphaEwald);
Peter Eastman's avatar
Peter Eastman committed
198
199
200
    RealOpenMM SQRT_PI                  = sqrt(PI_M);
    RealOpenMM TWO_PI                   = 2.0 * PI_M;
    RealOpenMM recipCoeff               = (RealOpenMM)(ONE_4PI_EPS0*4*PI_M/(periodicBoxSize[0] * periodicBoxSize[1] * periodicBoxSize[2]) /epsilon);
201

202
    RealOpenMM totalSelfEwaldEnergy     = 0.0;
203
204
    RealOpenMM realSpaceEwaldEnergy     = 0.0;
    RealOpenMM recipEnergy              = 0.0;
205
    RealOpenMM totalRecipEnergy         = 0.0;
206
    RealOpenMM vdwEnergy                = 0.0;
207
208
209
210
211

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

212
213
214
215
216
217
218
    if (includeReciprocal) {
        for( int atomID = 0; atomID < numberOfAtoms; atomID++ ){
            RealOpenMM selfEwaldEnergy       = (RealOpenMM) (ONE_4PI_EPS0*atomParameters[atomID][QIndex]*atomParameters[atomID][QIndex] * alphaEwald/SQRT_PI);
            totalSelfEwaldEnergy            -= selfEwaldEnergy;
            if( energyByAtom ){
                energyByAtom[atomID]        -= selfEwaldEnergy;
            }
219
        }
220
221
    }

222
223
224
225
    if( totalEnergy ){
        *totalEnergy += totalSelfEwaldEnergy;
    }

226
227
228
// **************************************************************************************
// RECIPROCAL SPACE EWALD ENERGY AND FORCES
// **************************************************************************************
229
    // PME
230

231
  if (pme && includeReciprocal) {
232
233
    pme_t          pmedata; /* abstract handle for PME data */
    RealOpenMM virial[3][3];
234

235
    pme_init(&pmedata,alphaEwald,numberOfAtoms,meshDim,5,1);
236

237
    pme_exec(pmedata,atomCoordinates,forces,atomParameters,periodicBoxSize,&recipEnergy,virial);
238

239
    if( totalEnergy )
240
241
242
243
244
245
       *totalEnergy += recipEnergy;

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

246
        pme_destroy(pmedata);
247
248
249
250
  }

    // Ewald method

251
  else if (ewald && includeReciprocal) {
252
253

    // setup reciprocal box
254
255
256
257

           RealOpenMM recipBoxSize[3] = { TWO_PI / periodicBoxSize[0], TWO_PI / periodicBoxSize[1], TWO_PI / periodicBoxSize[2]};


258
    // setup K-vectors
259

260
  #define EIR(x, y, z) eir[(x)*numberOfAtoms*3+(y)*3+z]
261
262
263
  vector<d_complex> eir(kmax*numberOfAtoms*3);
  vector<d_complex> tab_xy(numberOfAtoms);
  vector<d_complex> tab_qxyz(numberOfAtoms);
264
265
266
267
268
269
270

  if (kmax < 1) {
      std::stringstream message;
      message << " kmax < 1 , Aborting" << std::endl;
      SimTKOpenMMLog::printError( message );
  }

271
272
  for(int i = 0; (i < numberOfAtoms); i++) {
    for(int m = 0; (m < 3); m++)
273
      EIR(0, i, m) = d_complex(1,0);
274

275
    for(int m=0; (m<3); m++)
276
277
      EIR(1, i, m) = d_complex(cos(atomCoordinates[i][m]*recipBoxSize[m]),
                               sin(atomCoordinates[i][m]*recipBoxSize[m]));
278

279
280
    for(int j=2; (j<kmax); j++)
      for(int m=0; (m<3); m++)
281
        EIR(j, i, m) = EIR(j-1, i, m) * EIR(1, i, m);
282
283
  }

284
285
    // calculate reciprocal space energy and forces

286
287
288
289
290
291
292
293
294
295
296
297
    int lowry = 0;
    int lowrz = 1;

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

      RealOpenMM kx = rx * recipBoxSize[0];

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

        RealOpenMM ky = ry * recipBoxSize[1];

        if(ry >= 0) {
298
          for(int n = 0; n < numberOfAtoms; n++)
299
            tab_xy[n] = EIR(rx, n, 0) * EIR(ry, n, 1);
300
301
302
        }

        else {
303
          for(int n = 0; n < numberOfAtoms; n++)
304
            tab_xy[n]= EIR(rx, n, 0) * conj (EIR(-ry, n, 1));
305
306
307
308
309
        }

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

          if( rz >= 0) {
310
311
           for( int n = 0; n < numberOfAtoms; n++)
             tab_qxyz[n] = atomParameters[n][QIndex] * (tab_xy[n] * EIR(rz, n, 2));
312
313
314
          }

          else {
315
316
            for( int n = 0; n < numberOfAtoms; n++)
              tab_qxyz[n] = atomParameters[n][QIndex] * (tab_xy[n] * conj(EIR(-rz, n, 2)));
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
          }

          RealOpenMM cs = 0.0f;
          RealOpenMM ss = 0.0f;

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

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

          for(int n = 0; n < numberOfAtoms; n++) {
            RealOpenMM force = ak * (cs * tab_qxyz[n].imag() - ss * tab_qxyz[n].real());
333
334
335
            forces[n][0] += 2 * recipCoeff * force * kx ;
            forces[n][1] += 2 * recipCoeff * force * ky ;
            forces[n][2] += 2 * recipCoeff * force * kz ;
336
          }
337

338
339
          recipEnergy       = recipCoeff * ak * ( cs * cs + ss * ss);
          totalRecipEnergy += recipEnergy;
340
341
342
343
344
345
346
347

          if( totalEnergy )
             *totalEnergy += recipEnergy;

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

348
349
350
351
352
          lowrz = 1 - numRz;
        }
        lowry = 1 - numRy;
      }
    }
353
  }
354
355
356
357
358

// **************************************************************************************
// SHORT-RANGE ENERGY AND FORCES
// **************************************************************************************

359
360
    if (!includeDirect)
        return;
361
362
363
    RealOpenMM totalVdwEnergy            = 0.0f;
    RealOpenMM totalRealSpaceEwaldEnergy = 0.0f;

364
365
366
367
    for (int i = 0; i < (int) neighborList->size(); i++) {
       OpenMM::AtomPair pair = (*neighborList)[i];
       int ii = pair.first;
       int jj = pair.second;
368

369
370
       RealOpenMM deltaR[2][ReferenceForce::LastDeltaRIndex];
       ReferenceForce::getDeltaRPeriodic( atomCoordinates[jj], atomCoordinates[ii], periodicBoxSize, deltaR[0] );
371
372
       RealOpenMM r         = deltaR[0][ReferenceForce::RIndex];
       RealOpenMM inverseR  = one/(deltaR[0][ReferenceForce::RIndex]);
373
374
375
376
377
378
       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);
       }
379
       RealOpenMM alphaR    = alphaEwald * r;
380
381


Peter Eastman's avatar
Peter Eastman committed
382
383
       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 ));
384

385
386
387
388
389
       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];
390
391
392
393
394
395
                  dEdR     += switchValue*eps*( twelve*sig6 - six )*sig6*inverseR*inverseR;
       vdwEnergy = eps*(sig6-one)*sig6;
       if (useSwitch) {
           dEdR -= vdwEnergy*switchDeriv*inverseR;
           vdwEnergy *= switchValue;
       }
396

397
       // accumulate forces
398

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

405
       // accumulate energies
406

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

409
410
       totalVdwEnergy             += vdwEnergy;
       totalRealSpaceEwaldEnergy  += realSpaceEwaldEnergy;
411

412
413
414
415
        if( energyByAtom ){
           energyByAtom[ii] += realSpaceEwaldEnergy + vdwEnergy;
           energyByAtom[jj] += realSpaceEwaldEnergy + vdwEnergy;
        }
416

417
    }
418

419
420
421
    if( totalEnergy )
        *totalEnergy += totalRealSpaceEwaldEnergy + totalVdwEnergy;

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

424
    RealOpenMM totalExclusionEnergy = 0.0f;
Peter Eastman's avatar
Peter Eastman committed
425
426
427
428
429
430
431
    for (int i = 0; i < numberOfAtoms; i++)
        for (int j = 1; j <= exclusions[i][0]; j++)
            if (exclusions[i][j] > i) {
               int ii = i;
               int jj = exclusions[i][j];

               RealOpenMM deltaR[2][ReferenceForce::LastDeltaRIndex];
432
               ReferenceForce::getDeltaR( atomCoordinates[jj], atomCoordinates[ii], deltaR[0] );
Peter Eastman's avatar
Peter Eastman committed
433
434
435
               RealOpenMM r         = deltaR[0][ReferenceForce::RIndex];
               RealOpenMM inverseR  = one/(deltaR[0][ReferenceForce::RIndex]);
               RealOpenMM alphaR    = alphaEwald * r;
Peter Eastman's avatar
Peter Eastman committed
436
437
               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 ));
Peter Eastman's avatar
Peter Eastman committed
438
439
440
441
442
443
444
445
446
447
448

               // 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

449
               realSpaceEwaldEnergy = (RealOpenMM) (ONE_4PI_EPS0*atomParameters[ii][QIndex]*atomParameters[jj][QIndex]*inverseR*erf(alphaR));
Peter Eastman's avatar
Peter Eastman committed
450

451
452
               totalExclusionEnergy += realSpaceEwaldEnergy;
               if( energyByAtom ){
Peter Eastman's avatar
Peter Eastman committed
453
454
                   energyByAtom[ii] -= realSpaceEwaldEnergy;
                   energyByAtom[jj] -= realSpaceEwaldEnergy;
455
               }
Peter Eastman's avatar
Peter Eastman committed
456
457
            }

458
459
    if( totalEnergy )
        *totalEnergy -= totalExclusionEnergy;
460
461
462
}


463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
/**---------------------------------------------------------------------------------------

   Calculate LJ Coulomb pair ixn

   @param numberOfAtoms    number of atoms
   @param atomCoordinates  atom coordinates
   @param atomParameters   atom parameters                             atomParameters[atomIndex][paramterIndex]
   @param exclusions       atom exclusion indices                      exclusions[atomIndex][atomToExcludeIndex]
                           exclusions[atomIndex][0] = number of exclusions
                           exclusions[atomIndex][1-no.] = atom indices of atoms to excluded from
                           interacting w/ atom atomIndex
   @param fixedParameters  non atom parameters (not currently used)
   @param forces           force array (forces added)
   @param energyByAtom     atom energy
   @param totalEnergy      total energy
478
479
   @param includeDirect      true if direct space interactions should be included
   @param includeReciprocal  true if reciprocal space interactions should be included
480
481

   --------------------------------------------------------------------------------------- */
482

483
void ReferenceLJCoulombIxn::calculatePairIxn(int numberOfAtoms, vector<RealVec>& atomCoordinates,
484
                                             RealOpenMM** atomParameters, int** exclusions,
485
                                             RealOpenMM* fixedParameters, vector<RealVec>& forces,
486
                                             RealOpenMM* energyByAtom, RealOpenMM* totalEnergy, bool includeDirect, bool includeReciprocal) const {
487

488
489
490
491
492
493
494
   if (ewald || pme) {
       calculateEwaldIxn(numberOfAtoms, atomCoordinates, atomParameters, exclusions, fixedParameters, forces, energyByAtom,
               totalEnergy, includeDirect, includeReciprocal);
       return;
   }
   if (!includeDirect)
       return;
495
   if (cutoff) {
496
       for (int i = 0; i < (int) neighborList->size(); i++) {
497
498
499
500
501
502
           OpenMM::AtomPair pair = (*neighborList)[i];
           calculateOneIxn(pair.first, pair.second, atomCoordinates, atomParameters, forces, energyByAtom, totalEnergy);
       }
   }
   else {
       // allocate and initialize exclusion array
503

504
505
506
507
       int* exclusionIndices = new int[numberOfAtoms];
       for( int ii = 0; ii < numberOfAtoms; ii++ ){
          exclusionIndices[ii] = -1;
       }
508

509
       for( int ii = 0; ii < numberOfAtoms; ii++ ){
510

511
          // set exclusions
512

513
514
515
          for( int jj = 1; jj <= exclusions[ii][0]; jj++ ){
             exclusionIndices[exclusions[ii][jj]] = ii;
          }
516

517
          // loop over atom pairs
518

519
          for( int jj = ii+1; jj < numberOfAtoms; jj++ ){
520

521
522
523
524
525
             if( exclusionIndices[jj] != ii ){
                 calculateOneIxn(ii, jj, atomCoordinates, atomParameters, forces, energyByAtom, totalEnergy);
             }
          }
       }
526

527
528
529
       delete[] exclusionIndices;
   }
}
530

531
  /**---------------------------------------------------------------------------------------
532

533
     Calculate LJ Coulomb pair ixn between two atoms
534

535
536
537
538
539
540
541
     @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
542

543
     --------------------------------------------------------------------------------------- */
544

545
void ReferenceLJCoulombIxn::calculateOneIxn( int ii, int jj, vector<RealVec>& atomCoordinates,
546
                        RealOpenMM** atomParameters, vector<RealVec>& forces,
547
                        RealOpenMM* energyByAtom, RealOpenMM* totalEnergy ) const {
548

549
    // ---------------------------------------------------------------------------------------
550

551
    static const std::string methodName = "\nReferenceLJCoulombIxn::calculateOneIxn";
552

553
    // ---------------------------------------------------------------------------------------
554

555
    // constants -- reduce Visual Studio warnings regarding conversions between float & double
556

557
558
559
560
561
562
563
    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;
564

565
    static const int threeI             = 3;
566

567
568
569
570
571
572
573
    static const int LastAtomIndex      = 2;

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

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

    if (periodic)
574
        ReferenceForce::getDeltaRPeriodic( atomCoordinates[jj], atomCoordinates[ii], periodicBoxSize, deltaR[0] );
575
    else
576
        ReferenceForce::getDeltaR( atomCoordinates[jj], atomCoordinates[ii], deltaR[0] );
577
578
579

    RealOpenMM r2        = deltaR[0][ReferenceForce::R2Index];
    RealOpenMM inverseR  = one/(deltaR[0][ReferenceForce::RIndex]);
580
581
582
583
584
585
586
587
588
    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);
        }
    }
589
590
591
592
593
594
    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];
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
    RealOpenMM dEdR      = switchValue*eps*( twelve*sig6 - six )*sig6;
    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);
610

611
612
613
614
615
616
617
618
619
620
    // 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

621
622
623
624
625
    if( totalEnergy )
       *totalEnergy += energy;
    if( energyByAtom ){
       energyByAtom[ii] += energy;
       energyByAtom[jj] += energy;
626
627
    }
  }
628