TestReferenceEwald.cpp 18.2 KB
Newer Older
1
2
3
4
5
6
7
8
/* -------------------------------------------------------------------------- *
 *                                   OpenMM                                   *
 * -------------------------------------------------------------------------- *
 * This is part of the OpenMM molecular simulation toolkit originating from   *
 * Simbios, the NIH National Center for Physics-Based Simulation of           *
 * Biological Structures at Stanford, funded under the NIH Roadmap for        *
 * Medical Research, grant U54 GM072970. See https://simtk.org.               *
 *                                                                            *
9
 * Portions copyright (c) 2008-2014 Stanford University and the Authors.      *
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
 * Authors: Peter Eastman                                                     *
 * Contributors:                                                              *
 *                                                                            *
 * 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.                                     *
 * -------------------------------------------------------------------------- */

/**
 * This tests the Eewald summation method reference implementation of NonbondedForce.
 */

36
#include "openmm/internal/AssertionUtilities.h"
37
#include "openmm/Context.h"
38
#include "ReferencePlatform.h"
39
40
41
#include "openmm/NonbondedForce.h"
#include "openmm/System.h"
#include "openmm/VerletIntegrator.h"
42
#include "SimTKOpenMMRealType.h"
43
#include "sfmt/SFMT.h"
44
#include "openmm/HarmonicBondForce.h"
45
46
47
48
49
50
#include <iostream>
#include <vector>

using namespace OpenMM;
using namespace std;

51
52
ReferencePlatform platform;

Peter Eastman's avatar
Peter Eastman committed
53
54
const double EWALD_TOL = 1e-5;
const double PME_TOL = 5e-5;
55
56
57
58
59
60
61
62
63

void testEwaldExact() {

//    Use a NaCl crystal to compare the calculated and Madelung energies

    const int numParticles = 1000;
    const double cutoff = 1.0;
    const double boxSize = 2.82;

64
    System system;
65
    for (int i = 0; i < numParticles/2; i++)
66
        system.addParticle(22.99);
67
    for (int i = 0; i < numParticles/2; i++)
68
69
70
        system.addParticle(35.45);
    VerletIntegrator integrator(0.01);
    NonbondedForce* nonbonded = new NonbondedForce();
71
    for (int i = 0; i < numParticles/2; i++)
72
        nonbonded->addParticle(1.0, 1.0,0.0);
73
    for (int i = 0; i < numParticles/2; i++)
74
        nonbonded->addParticle(-1.0, 1.0,0.0);
75
    nonbonded->setNonbondedMethod(NonbondedForce::Ewald);
76
    nonbonded->setCutoffDistance(cutoff);
77
    system.setDefaultPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
Peter Eastman's avatar
Peter Eastman committed
78
    nonbonded->setEwaldErrorTolerance(EWALD_TOL);
79
80
    system.addForce(nonbonded);
    Context context(system, integrator, platform);
81
    vector<Vec3> positions(numParticles);
82
83
    #include "nacl_crystal.dat"
    context.setPositions(positions);
84

85
86
    State state = context.getState(State::Forces | State::Energy);
    const vector<Vec3>& forces = state.getForces();
87

88
//   The potential energy of an ion in a crystal is 
89
//   E = - (M*e^2/ 4*pi*epsilon0*a0),
90
//   where 
91
92
93
94
//   M            :    Madelung constant (dimensionless, for FCC cells such as NaCl it is 1.7476)
//   e            :    1.6022 × 10−19 C
//   4*pi*epsilon0:     1.112 × 10−10 C²/(J m)
//   a0           :    0.282 x 10-9 m (perfect cell)
95
96
97
// 
//   E is then the energy per pair of ions, so for our case
//   E has to be divided by 2 (per ion), multiplied by N(avogadro), multiplied by number of particles, and divided by 1000 for kJ
98
    double exactEnergy        = - (1.7476 * 1.6022e-19 * 1.6022e-19  * 6.02214e+23 * numParticles) / (1.112e-10 * 0.282e-9 * 2 * 1000);
99
100
    //cout << "exact\t\t: " << exactEnergy << endl;
    //cout << "calc\t\t: " << state.getPotentialEnergy() << endl;
Peter Eastman's avatar
Peter Eastman committed
101
    ASSERT_EQUAL_TOL(exactEnergy, state.getPotentialEnergy(), 100*EWALD_TOL);
102

103
104
}

105
106
void testEwaldPME() {

107
    double tol = 1e-5;
108
109
110
    const double boxSize = 3.00646;
    const double cutoff = 1.2;
    const int numParticles = 894;
111
112

//      Use amorphous NaCl system
113
//      The particles are simple charges, no VdW interactions
114
115
116
117
118
119
120
121
122
123
124
125
126
127

    System system;
    for (int i = 0; i < numParticles/2; i++)
        system.addParticle(22.99);
    for (int i = 0; i < numParticles/2; i++)
        system.addParticle(35.45);
    VerletIntegrator integrator(0.01);
    NonbondedForce* nonbonded = new NonbondedForce();
    for (int i = 0; i < numParticles/2; i++)
        nonbonded->addParticle(1.0, 1.0,0.0);
    for (int i = 0; i < numParticles/2; i++)
        nonbonded->addParticle(-1.0, 1.0,0.0);
    nonbonded->setNonbondedMethod(NonbondedForce::Ewald);
    nonbonded->setCutoffDistance(cutoff);
128
    system.setDefaultPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
Peter Eastman's avatar
Peter Eastman committed
129
    nonbonded->setEwaldErrorTolerance(EWALD_TOL);
130
131
132
133
134
135
136
137
138
    system.addForce(nonbonded);
    Context context(system, integrator, platform);
    vector<Vec3> positions(numParticles);
    #include "nacl_amorph.dat"
    context.setPositions(positions);

    State state1 = context.getState(State::Forces | State::Energy);
    const vector<Vec3>& forces1 = state1.getForces();

139
//   (1)   CHECK EXACT VALUE OF EWALD ENERGY (Against Gromacs output)
140

141
    tol = 1e-4;
142
    ASSERT_EQUAL_TOL(-3.82047e+05, state1.getPotentialEnergy(), tol);
143
144

//   (2)   CHECK WHETHER THE EWALD FORCES ARE THE SAME AS THE GROMACS OUTPUT
145
//         these are forces for alpha: 2.82756, kmax(x/y/z) = 11
146
  
147
148
    tol = 1e-2;
//    #include "nacl_amorph_GromacsForcesEwald.dat"
149
150
151

//   (3)   CHECK SELF-CONSISTENCY

152
    // Take a small step in the direction of the energy gradient.
153
154
155
156
157
158
159
160
161

    double norm = 0.0;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 f = state1.getForces()[i];
        norm += f[0]*f[0] + f[1]*f[1] + f[2]*f[2];
    }


    norm = std::sqrt(norm);
Peter Eastman's avatar
Peter Eastman committed
162
    const double delta = 1e-2;
163
164
165
166
167
168
169
170
171
172
173
    double step = delta/norm;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 p = positions[i];
        Vec3 f = state1.getForces()[i];
        positions[i] = Vec3(p[0]-f[0]*step, p[1]-f[1]*step, p[2]-f[2]*step);
    }
    context.setPositions(positions);
    
    // See whether the potential energy changed by the expected amount.
    
    State state2 = context.getState(State::Energy);
Peter Eastman's avatar
Peter Eastman committed
174
    ASSERT_EQUAL_TOL(norm, (state2.getPotentialEnergy()-state1.getPotentialEnergy())/delta, fabs(EWALD_TOL*state2.getPotentialEnergy()/(state2.getPotentialEnergy()-state1.getPotentialEnergy())))
175

176
//   (4)   CHECK EXACT VALUE OF PME ENERGY 
177
178

    nonbonded->setNonbondedMethod(NonbondedForce::PME);
Peter Eastman's avatar
Peter Eastman committed
179
    nonbonded->setEwaldErrorTolerance(PME_TOL);
180
181
182
183
    context.reinitialize();
    #include "nacl_amorph.dat"
    context.setPositions(positions);
    State state3 = context.getState(State::Forces | State::Energy);
184

185
//  Gromacs PME energy for the same mesh
186
    tol = 1e-4;
187
188
    ASSERT_EQUAL_TOL(-3.82047e+05, state3.getPotentialEnergy(), tol);

189
//   (5) CHECK WHETHER PME FORCES ARE THE SAME AS THE GROMACS OUTPUT USING EWALD
190

191
    tol = 1e-1;
192
//    #include "nacl_amorph_GromacsForcesEwald.dat"
193
194
195

//   (6) CHECK PME FOR SELF-CONSISTENCY

196
    // Take a small step in the direction of the energy gradient.
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
    
    norm = 0.0;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 f = state3.getForces()[i];
        norm += f[0]*f[0] + f[1]*f[1] + f[2]*f[2];
    }
    norm = std::sqrt(norm);
    step = delta/norm;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 p = positions[i];
        Vec3 f = state3.getForces()[i];
        positions[i] = Vec3(p[0]-f[0]*step, p[1]-f[1]*step, p[2]-f[2]*step);
    }
    context.setPositions(positions);
    
    // See whether the potential energy changed by the expected amount.
    
    State state4 = context.getState(State::Energy);
Peter Eastman's avatar
Peter Eastman committed
215
    ASSERT_EQUAL_TOL(norm, (state4.getPotentialEnergy()-state3.getPotentialEnergy())/delta, fabs(PME_TOL*state4.getPotentialEnergy()/(state4.getPotentialEnergy()-state3.getPotentialEnergy())))
216
217
}

218
219
220
221
222
223
224
225
226
227
228
void testEwald2Ions() {
    System system;
    system.addParticle(1.0);
    system.addParticle(1.0);
    VerletIntegrator integrator(0.01);
    NonbondedForce* nonbonded = new NonbondedForce();
    nonbonded->addParticle(1.0, 1, 0);
    nonbonded->addParticle(-1.0, 1, 0);
    nonbonded->setNonbondedMethod(NonbondedForce::Ewald);
    const double cutoff = 2.0;
    nonbonded->setCutoffDistance(cutoff);
229
    system.setDefaultPeriodicBoxVectors(Vec3(6, 0, 0), Vec3(0, 6, 0), Vec3(0, 0, 6));
Peter Eastman's avatar
Peter Eastman committed
230
    nonbonded->setEwaldErrorTolerance(EWALD_TOL);
231
232
233
234
235
236
237
238
    system.addForce(nonbonded);
    Context context(system, integrator, platform);
    vector<Vec3> positions(2);
    positions[0] = Vec3(3.048000,2.764000,3.156000);
    positions[1] = Vec3(2.809000,2.888000,2.571000);
    context.setPositions(positions);
    State state = context.getState(State::Forces | State::Energy);
    const vector<Vec3>& forces = state.getForces();
Peter Eastman's avatar
Peter Eastman committed
239
240
241
    ASSERT_EQUAL_VEC(Vec3(-123.711, 64.1877, -302.716), forces[0], 10*EWALD_TOL);
    ASSERT_EQUAL_VEC(Vec3(123.711, -64.1877, 302.716), forces[1], 10*EWALD_TOL);
    ASSERT_EQUAL_TOL(-217.276, state.getPotentialEnergy(), 10*EWALD_TOL);
242
243
244
245
}

void testWaterSystem() {
    System system;
246
247
248
    static int numParticles = 648;
    const double boxSize = 1.86206;

249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
    for (int i = 0 ; i < numParticles ; i++)
    {
       system.addParticle(1.0);
    }
    VerletIntegrator integrator(0.01);
    NonbondedForce* nonbonded = new NonbondedForce();
    for (int i = 0 ; i < numParticles/3 ; i++)
    {
      nonbonded->addParticle(-0.82, 1, 0);
      nonbonded->addParticle(0.41, 1, 0);
      nonbonded->addParticle(0.41, 1, 0);
    }
    nonbonded->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
    const double cutoff = 0.8;
    nonbonded->setCutoffDistance(cutoff);
264
    system.setDefaultPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
Peter Eastman's avatar
Peter Eastman committed
265
    nonbonded->setEwaldErrorTolerance(EWALD_TOL);
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
    system.addForce(nonbonded);
    Context context(system, integrator, platform);
    vector<Vec3> positions(numParticles);
    #include "water.dat"
    context.setPositions(positions);
    State state1 = context.getState(State::Forces | State::Energy);
    const vector<Vec3>& forces = state1.getForces();

// Take a small step in the direction of the energy gradient.
    
    double norm = 0.0;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 f = state1.getForces()[i];
        norm += f[0]*f[0] + f[1]*f[1] + f[2]*f[2];
    }


    norm = std::sqrt(norm);
    const double delta = 1e-3;
    double step = delta/norm;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 p = positions[i];
        Vec3 f = state1.getForces()[i];
        positions[i] = Vec3(p[0]-f[0]*step, p[1]-f[1]*step, p[2]-f[2]*step);
    }
    context.setPositions(positions);
    
    // See whether the potential energy changed by the expected amount.
    
    nonbonded->setNonbondedMethod(NonbondedForce::Ewald);
    State state2 = context.getState(State::Energy);
    ASSERT_EQUAL_TOL(norm, (state2.getPotentialEnergy()-state1.getPotentialEnergy())/delta, 0.01)


}
301

302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
void testTriclinic() {
    // Create a triclinic box containing eight particles.

    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(2.5, 0, 0), Vec3(0.5, 3.0, 0), Vec3(0.7, 0.9, 3.5));
    for (int i = 0; i < 8; i++)
        system.addParticle(1.0);
    NonbondedForce* force = new NonbondedForce();
    system.addForce(force);
    force->setNonbondedMethod(NonbondedForce::PME);
    force->setCutoffDistance(1.0);
    force->setPMEParameters(3.45891, 32, 40, 48);
    for (int i = 0; i < 4; i++)
        force->addParticle(-1, 0.440104, 0.4184); // Cl parameters
    for (int i = 0; i < 4; i++)
        force->addParticle(1, 0.332840, 0.0115897); // Na parameters
    vector<Vec3> positions(8);
    positions[0] = Vec3(1.744, 2.788, 3.162);
    positions[1] = Vec3(1.048, 0.762, 2.340);
    positions[2] = Vec3(2.489, 1.570, 2.817);
    positions[3] = Vec3(1.027, 1.893, 3.271);
    positions[4] = Vec3(0.937, 0.825, 0.009);
    positions[5] = Vec3(2.290, 1.887, 3.352);
    positions[6] = Vec3(1.266, 1.111, 2.894);
    positions[7] = Vec3(0.933, 1.862, 3.490);

    // Compute the forces and energy.

    VerletIntegrator integ(0.001);
    Context context(system, integ, platform);
    context.setPositions(positions);
    State state = context.getState(State::Forces | State::Energy);

    // Compare them to values computed by Gromacs.

    double expectedEnergy = -963.370;
    vector<Vec3> expectedForce(8);
    expectedForce[0] = Vec3(4.25253e+01, -1.23503e+02, 1.22139e+02);
    expectedForce[1] = Vec3(9.74752e+01, 1.68213e+02, 1.93169e+02);
    expectedForce[2] = Vec3(-1.50348e+02, 1.29165e+02, 3.70435e+02);
    expectedForce[3] = Vec3(9.18644e+02, -3.52571e+00, -1.34772e+03);
    expectedForce[4] = Vec3(-1.61193e+02, 9.01528e+01, -7.12904e+01);
    expectedForce[5] = Vec3(2.82630e+02, 2.78029e+01, -3.72864e+02);
    expectedForce[6] = Vec3(-1.47454e+02, -2.14448e+02, -3.55789e+02);
    expectedForce[7] = Vec3(-8.82195e+02, -7.39132e+01, 1.46202e+03);
    for (int i = 0; i < 8; i++) {
        ASSERT_EQUAL_VEC(expectedForce[i], state.getForces()[i], 1e-4);
    }
    ASSERT_EQUAL_TOL(expectedEnergy, state.getPotentialEnergy(), 1e-4);
}

353
354
355
356
357
358
void testErrorTolerance(NonbondedForce::NonbondedMethod method) {
    // Create a cloud of random point charges.

    const int numParticles = 51;
    const double boxWidth = 5.0;
    System system;
359
    system.setDefaultPeriodicBoxVectors(Vec3(boxWidth, 0, 0), Vec3(0, boxWidth, 0), Vec3(0, 0, boxWidth));
360
361
362
    NonbondedForce* force = new NonbondedForce();
    system.addForce(force);
    vector<Vec3> positions(numParticles);
363
364
365
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);

366
367
368
    for (int i = 0; i < numParticles; i++) {
        system.addParticle(1.0);
        force->addParticle(-1.0+i*2.0/(numParticles-1), 1.0, 0.0);
369
        positions[i] = Vec3(boxWidth*genrand_real2(sfmt), boxWidth*genrand_real2(sfmt), boxWidth*genrand_real2(sfmt));
370
371
372
373
374
    }
    force->setNonbondedMethod(method);

    // For various values of the cutoff and error tolerance, see if the actual error is reasonable.

375
    for (double cutoff = 1.0; cutoff < boxWidth/2; cutoff *= 1.2) {
376
377
378
379
380
        force->setCutoffDistance(cutoff);
        vector<Vec3> refForces;
        double norm = 0.0;
        for (double tol = 5e-5; tol < 1e-3; tol *= 2.0) {
            force->setEwaldErrorTolerance(tol);
Peter Eastman's avatar
Peter Eastman committed
381
            VerletIntegrator integrator(0.01);
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
            Context context(system, integrator, platform);
            context.setPositions(positions);
            State state = context.getState(State::Forces);
            if (refForces.size() == 0) {
                refForces = state.getForces();
                for (int i = 0; i < numParticles; i++)
                    norm += refForces[i].dot(refForces[i]);
                norm = sqrt(norm);
            }
            else {
                double diff = 0.0;
                for (int i = 0; i < numParticles; i++) {
                    Vec3 delta = refForces[i]-state.getForces()[i];
                    diff += delta.dot(delta);
                }
                diff = sqrt(diff)/norm;
398
                ASSERT(diff < 2*tol);
399
400
401
402
            }
        }
    }
}
403

404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
void testPMEParameters() {
    // Create a cloud of random point charges.

    const int numParticles = 51;
    const double boxWidth = 4.7;
    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(boxWidth, 0, 0), Vec3(0, boxWidth, 0), Vec3(0, 0, boxWidth));
    NonbondedForce* force = new NonbondedForce();
    system.addForce(force);
    vector<Vec3> positions(numParticles);
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);

    for (int i = 0; i < numParticles; i++) {
        system.addParticle(1.0);
        force->addParticle(-1.0+i*2.0/(numParticles-1), 1.0, 0.0);
        positions[i] = Vec3(boxWidth*genrand_real2(sfmt), boxWidth*genrand_real2(sfmt), boxWidth*genrand_real2(sfmt));
    }
    force->setNonbondedMethod(NonbondedForce::PME);
    
    // Compute the energy with an error tolerance of 1e-3.

    force->setEwaldErrorTolerance(1e-3);
    VerletIntegrator integrator1(0.01);
    Context context1(system, integrator1, platform);
    context1.setPositions(positions);
    double energy1 = context1.getState(State::Energy).getPotentialEnergy();
    
    // Try again with an error tolerance of 1e-4.

    force->setEwaldErrorTolerance(1e-4);
    VerletIntegrator integrator2(0.01);
    Context context2(system, integrator2, platform);
    context2.setPositions(positions);
    double energy2 = context2.getState(State::Energy).getPotentialEnergy();
    
    // Now explicitly set the parameters.  These should match the values that were
    // used for tolerance 1e-3.

    force->setPMEParameters(2.49291157051793, 32, 32, 32);
    VerletIntegrator integrator3(0.01);
    Context context3(system, integrator3, platform);
    context3.setPositions(positions);
    double energy3 = context3.getState(State::Energy).getPotentialEnergy();
    ASSERT_EQUAL_TOL(energy1, energy3, 1e-6);
    ASSERT(fabs((energy1-energy2)/energy1) > 1e-5);
}

452
453
int main() {
    try {
454
     testEwaldExact();
455
     testEwaldPME();
456
457
//     testEwald2Ions();
//     testWaterSystem();
458
     testTriclinic();
459
460
     testErrorTolerance(NonbondedForce::Ewald);
     testErrorTolerance(NonbondedForce::PME);
461
     testPMEParameters();
462
463
464
465
466
467
468
469
    }
    catch(const exception& e) {
        cout << "exception: " << e.what() << endl;
        return 1;
    }
    cout << "Done" << endl;
    return 0;
}