TestReferenceEwald.cpp 14 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-2009 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
36
 * 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.
 */

#include "../../../tests/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 "../src/SimTKUtilities/SimTKOpenMMRealType.h"
43
#include "../src/sfmt/SFMT.h"
44
#include "openmm/HarmonicBondForce.h"
45
46
47
48
49
50
51
52
#include <iostream>
#include <vector>

using namespace OpenMM;
using namespace std;

const double TOL = 1e-5;

53
54
55
56
57
58
59
60
61

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;

62
63
    ReferencePlatform platform;
    System system;
64
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.setPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
78
79
80
    nonbonded->setEwaldErrorTolerance(TOL);
    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
91
92
93
94
//   where 
//   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
98
// 
//   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
   	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;
101
    ASSERT_EQUAL_TOL(exactEnergy, state.getPotentialEnergy(), 100*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
128

    ReferencePlatform platform;
    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);
129
    system.setPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
130
131
132
133
134
135
136
137
138
139
    nonbonded->setEwaldErrorTolerance(TOL);
    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();

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

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

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

//   (3)   CHECK SELF-CONSISTENCY

153
    // Take a small step in the direction of the energy gradient.
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173

    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.
    
174
    tol = 1e-2;
175
    State state2 = context.getState(State::Energy);
176
    ASSERT_EQUAL_TOL(norm, (state2.getPotentialEnergy()-state1.getPotentialEnergy())/delta, tol)
177

178
//   (4)   CHECK EXACT VALUE OF PME ENERGY 
179
180
181
182
183
184

    nonbonded->setNonbondedMethod(NonbondedForce::PME);
    context.reinitialize();
    #include "nacl_amorph.dat"
    context.setPositions(positions);
    State state3 = context.getState(State::Forces | State::Energy);
185

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

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

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

//   (6) CHECK PME FOR SELF-CONSISTENCY

197
    // Take a small step in the direction of the energy gradient.
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
    
    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);
216
    tol = 1e-2;
217
    ASSERT_EQUAL_TOL(norm, (state4.getPotentialEnergy()-state3.getPotentialEnergy())/delta, tol)
218
219
220

}

221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
void testEwald2Ions() {
    ReferencePlatform platform;
    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);
    system.setPeriodicBoxVectors(Vec3(6, 0, 0), Vec3(0, 6, 0), Vec3(0, 0, 6));
    nonbonded->setEwaldErrorTolerance(TOL);
    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();
    ASSERT_EQUAL_VEC(Vec3(-123.711, 64.1877, -302.716), forces[0], 10*TOL);
    ASSERT_EQUAL_VEC(Vec3(123.711, -64.1877, 302.716), forces[1], 10*TOL);
    ASSERT_EQUAL_TOL(-217.276, state.getPotentialEnergy(), 10*TOL);
}

void testWaterSystem() {
    ReferencePlatform platform;
    System system;
251
252
253
    static int numParticles = 648;
    const double boxSize = 1.86206;

254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
    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);
269
    system.setPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
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
301
302
303
304
305
    nonbonded->setEwaldErrorTolerance(TOL);
    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)


}
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
353
354
355
void testErrorTolerance(NonbondedForce::NonbondedMethod method) {
    // Create a cloud of random point charges.

    const int numParticles = 51;
    const double boxWidth = 5.0;
    System system;
    system.setPeriodicBoxVectors(Vec3(boxWidth, 0, 0), Vec3(0, boxWidth, 0), Vec3(0, 0, boxWidth));
    NonbondedForce* force = new NonbondedForce();
    system.addForce(force);
    vector<Vec3> positions(numParticles);
    init_gen_rand(0);
    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(), boxWidth*genrand_real2(), boxWidth*genrand_real2());
    }
    force->setNonbondedMethod(method);
    ReferencePlatform platform;
    VerletIntegrator integrator(0.01);

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

    for (double cutoff = 1.0; cutoff < boxWidth/2; cutoff += 0.2) {
        force->setCutoffDistance(cutoff);
        vector<Vec3> refForces;
        double norm = 0.0;
        for (double tol = 5e-5; tol < 1e-3; tol *= 2.0) {
            force->setEwaldErrorTolerance(tol);
            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;
                ASSERT(diff < 5*tol);
            }
        }
    }
}
356
357
358

int main() {
    try {
359
     testEwaldExact();
360
     testEwaldPME();
361
362
//     testEwald2Ions();
//     testWaterSystem();
363
364
     testErrorTolerance(NonbondedForce::Ewald);
     testErrorTolerance(NonbondedForce::PME);
365
366
367
368
369
370
371
372
    }
    catch(const exception& e) {
        cout << "exception: " << e.what() << endl;
        return 1;
    }
    cout << "Done" << endl;
    return 0;
}