TestReferenceLangevinIntegrator.cpp 10.5 KB
Newer Older
1
2
3
4
5
6
7
8
9
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
/* -------------------------------------------------------------------------- *
 *                                   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.               *
 *                                                                            *
 * Portions copyright (c) 2008 Stanford University and the Authors.           *
 * 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 LangevinIntegrator.
 */

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

using namespace OpenMM;
using namespace std;

51
52
ReferencePlatform platform;

53
54
55
const double TOL = 1e-5;

void testSingleBond() {
56
57
58
    System system;
    system.addParticle(2.0);
    system.addParticle(2.0);
59
    LangevinIntegrator integrator(0, 0.1, 0.01);
60
61
    HarmonicBondForce* forceField = new HarmonicBondForce();
    forceField->addBond(0, 1, 1.5, 1);
62
    system.addForce(forceField);
63
    Context context(system, integrator, platform);
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
    vector<Vec3> positions(2);
    positions[0] = Vec3(-1, 0, 0);
    positions[1] = Vec3(1, 0, 0);
    context.setPositions(positions);
    
    // This is simply a damped harmonic oscillator, so compare it to the analytical solution.
    
    double freq = std::sqrt(1-0.05*0.05);
    for (int i = 0; i < 1000; ++i) {
        State state = context.getState(State::Positions | State::Velocities);
        double time = state.getTime();
        double expectedDist = 1.5+0.5*std::exp(-0.05*time)*std::cos(freq*time);
        ASSERT_EQUAL_VEC(Vec3(-0.5*expectedDist, 0, 0), state.getPositions()[0], 0.02);
        ASSERT_EQUAL_VEC(Vec3(0.5*expectedDist, 0, 0), state.getPositions()[1], 0.02);
        double expectedSpeed = -0.5*std::exp(-0.05*time)*(0.05*std::cos(freq*time)+freq*std::sin(freq*time));
        ASSERT_EQUAL_VEC(Vec3(-0.5*expectedSpeed, 0, 0), state.getVelocities()[0], 0.02);
        ASSERT_EQUAL_VEC(Vec3(0.5*expectedSpeed, 0, 0), state.getVelocities()[1], 0.02);
        integrator.step(1);
    }
    
    // Not set the friction to a tiny value and see if it conserves energy.
    
    integrator.setFriction(5e-5);
    context.setPositions(positions);
    State state = context.getState(State::Energy);
    double initialEnergy = state.getKineticEnergy()+state.getPotentialEnergy();
    for (int i = 0; i < 1000; ++i) {
        state = context.getState(State::Energy);
        double energy = state.getKineticEnergy()+state.getPotentialEnergy();
        ASSERT_EQUAL_TOL(initialEnergy, energy, 0.01);
        integrator.step(1);
    }
}

98
void testTemperature() {
Peter Eastman's avatar
Peter Eastman committed
99
    const int numParticles = 8;
100
    const double temp = 100.0;
101
    System system;
102
    LangevinIntegrator integrator(temp, 2.0, 0.01);
103
    NonbondedForce* forceField = new NonbondedForce();
Peter Eastman's avatar
Peter Eastman committed
104
    for (int i = 0; i < numParticles; ++i) {
105
        system.addParticle(2.0);
106
        forceField->addParticle((i%2 == 0 ? 1.0 : -1.0), 1.0, 5.0);
107
108
    }
    system.addForce(forceField);
109
    Context context(system, integrator, platform);
Peter Eastman's avatar
Peter Eastman committed
110
111
    vector<Vec3> positions(numParticles);
    for (int i = 0; i < numParticles; ++i)
112
113
114
115
116
117
118
119
120
121
        positions[i] = Vec3((i%2 == 0 ? 2 : -2), (i%4 < 2 ? 2 : -2), (i < 4 ? 2 : -2));
    context.setPositions(positions);
    
    // Let it equilibrate.
    
    integrator.step(10000);
    
    // Now run it for a while and see if the temperature is correct.
    
    double ke = 0.0;
122
    for (int i = 0; i < 10000; ++i) {
123
124
125
126
        State state = context.getState(State::Energy);
        ke += state.getKineticEnergy();
        integrator.step(1);
    }
127
    ke /= 10000;
Peter Eastman's avatar
Peter Eastman committed
128
    double expected = 0.5*numParticles*3*BOLTZ*temp;
129
    ASSERT_USUALLY_EQUAL_TOL(expected, ke, 6/std::sqrt(10000.0));
130
131
}

132
void testConstraints() {
Peter Eastman's avatar
Peter Eastman committed
133
    const int numParticles = 8;
134
    const double temp = 100.0;
135
    System system;
136
    LangevinIntegrator integrator(temp, 2.0, 0.01);
137
    integrator.setConstraintTolerance(1e-5);
138
    NonbondedForce* forceField = new NonbondedForce();
Peter Eastman's avatar
Peter Eastman committed
139
    for (int i = 0; i < numParticles; ++i) {
140
        system.addParticle(10.0);
141
        forceField->addParticle((i%2 == 0 ? 0.2 : -0.2), 0.5, 5.0);
142
    }
Peter Eastman's avatar
Peter Eastman committed
143
    for (int i = 0; i < numParticles-1; ++i)
144
        system.addConstraint(i, i+1, 1.0);
145
    system.addForce(forceField);
146
    Context context(system, integrator, platform);
Peter Eastman's avatar
Peter Eastman committed
147
148
    vector<Vec3> positions(numParticles);
    vector<Vec3> velocities(numParticles);
149
150
151
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);

Peter Eastman's avatar
Peter Eastman committed
152
    for (int i = 0; i < numParticles; ++i) {
153
        positions[i] = Vec3(i/2, (i+1)/2, 0);
154
        velocities[i] = Vec3(genrand_real2(sfmt)-0.5, genrand_real2(sfmt)-0.5, genrand_real2(sfmt)-0.5);
155
156
157
158
159
160
161
162
    }
    context.setPositions(positions);
    context.setVelocities(velocities);
    
    // Simulate it and see whether the constraints remain satisfied.
    
    for (int i = 0; i < 1000; ++i) {
        State state = context.getState(State::Positions);
Peter Eastman's avatar
Peter Eastman committed
163
        for (int j = 0; j < numParticles-1; ++j) {
164
165
            Vec3 p1 = state.getPositions()[j];
            Vec3 p2 = state.getPositions()[j+1];
166
            double dist = std::sqrt((p1[0]-p2[0])*(p1[0]-p2[0])+(p1[1]-p2[1])*(p1[1]-p2[1])+(p1[2]-p2[2])*(p1[2]-p2[2]));
167
            ASSERT_EQUAL_TOL(1.0, dist, 2e-5);
168
169
170
171
172
        }
        integrator.step(1);
    }
}

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
void testConstrainedMasslessParticles() {
    System system;
    system.addParticle(0.0);
    system.addParticle(1.0);
    system.addConstraint(0, 1, 1.5);
    vector<Vec3> positions(2);
    positions[0] = Vec3(-1, 0, 0);
    positions[1] = Vec3(1, 0, 0);
    LangevinIntegrator integrator(300.0, 2.0, 0.01);
    bool failed = false;
    try {
        // This should throw an exception.
        
        Context context(system, integrator, platform);
    }
    catch (exception& ex) {
        failed = true;
    }
    ASSERT(failed);
    
    // Now make both particles massless, which should work.
    
    system.setParticleMass(1, 0.0);
    Context context(system, integrator, platform);
    context.setPositions(positions);
    context.setVelocitiesToTemperature(300.0);
    integrator.step(1);
    State state = context.getState(State::Velocities | State::Positions);
    ASSERT_EQUAL(0.0, state.getVelocities()[0][0]);
}

204
205
206
207
void testRandomSeed() {
    const int numParticles = 8;
    const double temp = 100.0;
    const double collisionFreq = 10.0;
208
    System system;
209
    LangevinIntegrator integrator(temp, 2.0, 0.01);
210
    NonbondedForce* forceField = new NonbondedForce();
211
    for (int i = 0; i < numParticles; ++i) {
212
        system.addParticle(2.0);
213
        forceField->addParticle((i%2 == 0 ? 1.0 : -1.0), 1.0, 5.0);
214
215
216
217
218
219
220
221
222
223
224
225
    }
    system.addForce(forceField);
    vector<Vec3> positions(numParticles);
    vector<Vec3> velocities(numParticles);
    for (int i = 0; i < numParticles; ++i) {
        positions[i] = Vec3((i%2 == 0 ? 2 : -2), (i%4 < 2 ? 2 : -2), (i < 4 ? 2 : -2));
        velocities[i] = Vec3(0, 0, 0);
    }

    // Try twice with the same random seed.

    integrator.setRandomNumberSeed(5);
226
    Context context(system, integrator, platform);
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
    context.setPositions(positions);
    context.setVelocities(velocities);
    integrator.step(10);
    State state1 = context.getState(State::Positions);
    context.reinitialize();
    context.setPositions(positions);
    context.setVelocities(velocities);
    integrator.step(10);
    State state2 = context.getState(State::Positions);

    // Try twice with a different random seed.

    integrator.setRandomNumberSeed(10);
    context.reinitialize();
    context.setPositions(positions);
    context.setVelocities(velocities);
    integrator.step(10);
    State state3 = context.getState(State::Positions);
    context.reinitialize();
    context.setPositions(positions);
    context.setVelocities(velocities);
    integrator.step(10);
    State state4 = context.getState(State::Positions);

    // Compare the results.

    for (int i = 0; i < numParticles; i++) {
        for (int j = 0; j < 3; j++) {
            ASSERT(state1.getPositions()[i][j] == state2.getPositions()[i][j]);
            ASSERT(state3.getPositions()[i][j] == state4.getPositions()[i][j]);
            ASSERT(state1.getPositions()[i][j] != state3.getPositions()[i][j]);
        }
    }
}

262
263
264
int main() {
    try {
        testSingleBond();
265
        testTemperature();
266
        testConstraints();
267
        testConstrainedMasslessParticles();
268
        testRandomSeed();
269
270
271
272
273
274
275
276
    }
    catch(const exception& e) {
        cout << "exception: " << e.what() << endl;
        return 1;
    }
    cout << "Done" << endl;
    return 0;
}