TestReferenceRpmd.cpp 13.4 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) 2011-2013 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 Reference implementation of RPMDIntegrator.
 */

36
#include "openmm/internal/AssertionUtilities.h"
37
#include "openmm/CMMotionRemover.h"
38
39
#include "openmm/Context.h"
#include "openmm/HarmonicBondForce.h"
Peter Eastman's avatar
Peter Eastman committed
40
#include "openmm/NonbondedForce.h"
41
42
43
#include "openmm/Platform.h"
#include "openmm/System.h"
#include "openmm/RPMDIntegrator.h"
Peter Eastman's avatar
Peter Eastman committed
44
#include "openmm/VirtualSite.h"
45
#include "SimTKUtilities/SimTKOpenMMUtilities.h"
46
#include "sfmt/SFMT.h"
47
48
49
50
51
52
#include <iostream>
#include <vector>

using namespace OpenMM;
using namespace std;

Peter Eastman's avatar
Peter Eastman committed
53
54
void testFreeParticles() {
    const int numParticles = 100;
55
    const int numCopies = 30;
56
    const double temperature = 300.0;
Peter Eastman's avatar
Peter Eastman committed
57
    const double mass = 1.0;
58
    System system;
59
    for (int i = 0; i < numParticles; i++)
Peter Eastman's avatar
Peter Eastman committed
60
61
        system.addParticle(mass);
    RPMDIntegrator integ(numCopies, temperature, 10.0, 0.001);
Peter Eastman's avatar
Peter Eastman committed
62
63
    Platform& platform = Platform::getPlatformByName("Reference");
    Context context(system, integ, platform);
64
65
66
67
68
69
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    vector<Vec3> positions(numParticles);
    for (int i = 0; i < numCopies; i++)
    {
        for (int j = 0; j < numParticles; j++)
Peter Eastman's avatar
Peter Eastman committed
70
            positions[j] = Vec3(0.02*genrand_real2(sfmt), 0.02*genrand_real2(sfmt), 0.02*genrand_real2(sfmt));
71
72
        integ.setPositions(i, positions);
    }
Peter Eastman's avatar
Peter Eastman committed
73
74
    const int numSteps = 1000;
    integ.step(1000);
75
76
77
78
79
80
81
    vector<double> ke(numCopies, 0.0);
    vector<double> rg(numParticles, 0.0);
    const RealOpenMM hbar = 1.054571628e-34*AVOGADRO/(1000*1e-12);
    for (int i = 0; i < numSteps; i++) {
        integ.step(1);
        vector<State> state(numCopies);
        for (int j = 0; j < numCopies; j++)
82
            state[j] = integ.getState(j, State::Positions | State::Velocities, true);
83
84
        for (int j = 0; j < numParticles; j++) {
            double rg2 = 0.0;
Peter Eastman's avatar
Peter Eastman committed
85
86
87
            for (int k = 0; k < numCopies; k++) {
                Vec3 v = state[k].getVelocities()[j];
                ke[k] += 0.5*mass*v.dot(v);
88
89
90
91
                for (int m = 0; m < numCopies; m++) {
                    Vec3 delta = state[k].getPositions()[j]-state[m].getPositions()[j];
                    rg2 += delta.dot(delta);
                }
Peter Eastman's avatar
Peter Eastman committed
92
93
            }
            rg[j] += rg2/(2*numCopies*numCopies);
94
95
        }
    }
Peter Eastman's avatar
Peter Eastman committed
96
97
98
99
100
101
102
103
104
105
106
107
    double meanKE = 0.0;
    for (int i = 0; i < numCopies; i++)
        meanKE += ke[i];
    meanKE /= numSteps*numCopies;
    double expectedKE = 0.5*numCopies*numParticles*3*BOLTZ*temperature;
    ASSERT_USUALLY_EQUAL_TOL(expectedKE, meanKE, 1e-2);
    double meanRg2 = 0.0;
    for (int i = 0; i < numParticles; i++)
        meanRg2 += rg[i];
    meanRg2 /= numSteps*numParticles;
    double expectedRg = hbar/(2*sqrt(mass*BOLTZ*temperature));
    ASSERT_USUALLY_EQUAL_TOL(expectedRg, sqrt(meanRg2), 1e-3);
108
109
}

110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
Vec3 calcCM(const vector<Vec3>& values, System& system) {
    Vec3 cm;
    for (int j = 0; j < system.getNumParticles(); ++j) {
        cm[0] += values[j][0]*system.getParticleMass(j);
        cm[1] += values[j][1]*system.getParticleMass(j);
        cm[2] += values[j][2]*system.getParticleMass(j);
    }
    return cm;
}

void testCMMotionRemoval() {
    const int numParticles = 100;
    const int numCopies = 30;
    const double temperature = 300.0;
    const double mass = 1.0;
    System system;
    for (int i = 0; i < numParticles; i++)
        system.addParticle(mass);
    system.addForce(new CMMotionRemover());
    RPMDIntegrator integ(numCopies, temperature, 10.0, 0.001);
    Platform& platform = Platform::getPlatformByName("Reference");
    Context context(system, integ, platform);
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    vector<Vec3> positions(numParticles);
    for (int i = 0; i < numCopies; i++)
    {
        for (int j = 0; j < numParticles; j++)
            positions[j] = Vec3(0.02*genrand_real2(sfmt), 0.02*genrand_real2(sfmt), 0.02*genrand_real2(sfmt));
        Vec3 cmPos = calcCM(positions, system);
        for (int j = 0; j < numParticles; j++)
            positions[j] -= cmPos*(1/(mass*numParticles));
        integ.setPositions(i, positions);
    }
    
    // Make sure the CMMotionRemover is getting applied.
    
    for (int i = 0; i < 200; ++i) {
        integ.step(1);
        Vec3 pos;
        for (int j = 0; j < numCopies; j++) {
            State state = integ.getState(0, State::Positions | State::Velocities);
            pos += calcCM(state.getPositions(), system);
        }
        pos *= 1.0/numCopies;
        ASSERT_EQUAL_VEC(Vec3(), pos, 0.5);
    }
}

Peter Eastman's avatar
Peter Eastman committed
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
void testVirtualSites() {
    const int gridSize = 3;
    const int numMolecules = gridSize*gridSize*gridSize;
    const int numParticles = numMolecules*3;
    const int numCopies = 10;
    const double spacing = 2.0;
    const double cutoff = 3.0;
    const double boxSize = spacing*(gridSize+1);
    const double temperature = 300.0;
    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
    HarmonicBondForce* bonds = new HarmonicBondForce();
    system.addForce(bonds);
    NonbondedForce* nonbonded = new NonbondedForce();
    nonbonded->setCutoffDistance(cutoff);
    nonbonded->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
    system.addForce(nonbonded);

    // Create a cloud of molecules.

    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    vector<Vec3> positions(numParticles);
    for (int i = 0; i < numMolecules; i++) {
        system.addParticle(1.0);
        system.addParticle(1.0);
        system.addParticle(0.0);
        nonbonded->addParticle(-0.2, 0.2, 0.2);
        nonbonded->addParticle(0.1, 0.2, 0.2);
        nonbonded->addParticle(0.1, 0.2, 0.2);
        nonbonded->addException(3*i, 3*i+1, 0, 1, 0);
        nonbonded->addException(3*i, 3*i+2, 0, 1, 0);
        nonbonded->addException(3*i+1, 3*i+2, 0, 1, 0);
        bonds->addBond(3*i, 3*i+1, 1.0, 10000.0);
        system.setVirtualSite(3*i+2, new TwoParticleAverageSite(3*i, 3*i+1, 0.5, 0.5));
    }
    RPMDIntegrator integ(numCopies, temperature, 10.0, 0.001);
    Platform& platform = Platform::getPlatformByName("Reference");
    Context context(system, integ, platform);
    for (int copy = 0; copy < numCopies; copy++) {
        for (int i = 0; i < gridSize; i++)
            for (int j = 0; j < gridSize; j++)
                for (int k = 0; k < gridSize; k++) {
                    Vec3 pos = Vec3(spacing*(i+0.02*genrand_real2(sfmt)), spacing*(j+0.02*genrand_real2(sfmt)), spacing*(k+0.02*genrand_real2(sfmt)));
                    int index = k+gridSize*(j+gridSize*i);
                    positions[3*index] = pos;
                    positions[3*index+1] = Vec3(pos[0]+1.0, pos[1], pos[2]);
                    positions[3*index+2] = Vec3();
                }
        integ.setPositions(copy, positions);
    }

    // Check the temperature and virtual site locations.
    
    const int numSteps = 1000;
    integ.step(1000);
    vector<double> ke(numCopies, 0.0);
    for (int i = 0; i < numSteps; i++) {
        integ.step(1);
        vector<State> state(numCopies);
        for (int j = 0; j < numCopies; j++) {
            state[j] = integ.getState(j, State::Positions | State::Velocities | State::Forces, true);
            const vector<Vec3>& pos = state[j].getPositions();
            for (int k = 0; k < numMolecules; k++)
                ASSERT_EQUAL_VEC((pos[3*k]+pos[3*k+1])*0.5, pos[3*k+2], 1e-5);
        }
        for (int j = 0; j < numParticles; j++) {
            for (int k = 0; k < numCopies; k++) {
                Vec3 v = state[k].getVelocities()[j];
                ke[k] += 0.5*system.getParticleMass(j)*v.dot(v);
            }
        }
    }
    double meanKE = 0.0;
    for (int i = 0; i < numCopies; i++)
        meanKE += ke[i];
    meanKE /= numSteps*numCopies;
    double expectedKE = 0.5*numCopies*(2*numMolecules)*3*BOLTZ*temperature;
    ASSERT_USUALLY_EQUAL_TOL(expectedKE, meanKE, 1e-2);
}

240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
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
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
void testContractions() {
    const int gridSize = 3;
    const int numMolecules = gridSize*gridSize*gridSize;
    const int numParticles = numMolecules*2;
    const int numCopies = 10;
    const double spacing = 2.0;
    const double cutoff = 3.0;
    const double boxSize = spacing*(gridSize+1);
    const double temperature = 300.0;
    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
    HarmonicBondForce* bonds = new HarmonicBondForce();
    system.addForce(bonds);
    NonbondedForce* nonbonded = new NonbondedForce();
    nonbonded->setCutoffDistance(cutoff);
    nonbonded->setNonbondedMethod(NonbondedForce::PME);
    nonbonded->setForceGroup(1);
    nonbonded->setReciprocalSpaceForceGroup(2);
    system.addForce(nonbonded);

    // Create a cloud of molecules.

    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    vector<Vec3> positions(numParticles);
    for (int i = 0; i < numMolecules; i++) {
        system.addParticle(1.0);
        system.addParticle(1.0);
        nonbonded->addParticle(-0.2, 0.2, 0.2);
        nonbonded->addParticle(0.2, 0.2, 0.2);
        nonbonded->addException(2*i, 2*i+1, 0, 1, 0);
        bonds->addBond(2*i, 2*i+1, 1.0, 10000.0);
    }
    map<int, int> contractions;
    contractions[1] = 3;
    contractions[2] = 1;
    RPMDIntegrator integ(numCopies, temperature, 10.0, 0.001, contractions);
    Platform& platform = Platform::getPlatformByName("Reference");
    Context context(system, integ, platform);
    for (int copy = 0; copy < numCopies; copy++) {
        for (int i = 0; i < gridSize; i++)
            for (int j = 0; j < gridSize; j++)
                for (int k = 0; k < gridSize; k++) {
                    Vec3 pos = Vec3(spacing*(i+0.02*genrand_real2(sfmt)), spacing*(j+0.02*genrand_real2(sfmt)), spacing*(k+0.02*genrand_real2(sfmt)));
                    int index = k+gridSize*(j+gridSize*i);
                    positions[2*index] = pos;
                    positions[2*index+1] = Vec3(pos[0]+1.0, pos[1], pos[2]);
                }
        integ.setPositions(copy, positions);
    }

    // Check the temperature.
    
    const int numSteps = 1000;
    integ.step(1000);
    vector<double> ke(numCopies, 0.0);
    for (int i = 0; i < numSteps; i++) {
        integ.step(1);
        vector<State> state(numCopies);
        for (int j = 0; j < numCopies; j++)
            state[j] = integ.getState(j, State::Velocities, true);
        for (int j = 0; j < numParticles; j++) {
            for (int k = 0; k < numCopies; k++) {
                Vec3 v = state[k].getVelocities()[j];
                ke[k] += 0.5*system.getParticleMass(j)*v.dot(v);
            }
        }
    }
    double meanKE = 0.0;
    for (int i = 0; i < numCopies; i++)
        meanKE += ke[i];
    meanKE /= numSteps*numCopies;
    double expectedKE = 0.5*numCopies*numParticles*3*BOLTZ*temperature;
    ASSERT_USUALLY_EQUAL_TOL(expectedKE, meanKE, 1e-2);
}

316
317
int main() {
    try {
Peter Eastman's avatar
Peter Eastman committed
318
        testFreeParticles();
319
        testCMMotionRemoval();
Peter Eastman's avatar
Peter Eastman committed
320
        testVirtualSites();
321
        testContractions();
322
323
324
325
326
327
328
329
330
    }
    catch(const std::exception& e) {
        std::cout << "exception: " << e.what() << std::endl;
        std::cout << "FAIL - ERROR.  Test failed." << std::endl;
        return 1;
    }
    std::cout << "Done" << std::endl;
    return 0;
}