TestOpenCLCheckpoints.cpp 10 KB
Newer Older
Peter Eastman's avatar
Peter Eastman committed
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.               *
 *                                                                            *
Peter Eastman's avatar
Peter Eastman committed
9
 * Portions copyright (c) 2012-2013 Stanford University and the Authors.      *
Peter Eastman's avatar
Peter Eastman committed
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
37
38
39
40
41
42
43
44
45
46
47
48
49
50
 * 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 creating and loading checkpoints with the OpenCL platform.
 */

#include "OpenCLPlatform.h"
#include "openmm/internal/AssertionUtilities.h"
#include "openmm/AndersenThermostat.h"
#include "openmm/Context.h"
#include "openmm/NonbondedForce.h"
#include "openmm/System.h"
#include "openmm/VerletIntegrator.h"
#include "sfmt/SFMT.h"
#include <iostream>
#include <sstream>
#include <vector>

using namespace OpenMM;
using namespace std;

51
static OpenCLPlatform platform;
52

Peter Eastman's avatar
Peter Eastman committed
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
const double TOL = 1e-5;

void compareStates(State& s1, State& s2) {
    ASSERT_EQUAL_TOL(s1.getTime(), s2.getTime(), TOL);
    int numParticles = s1.getPositions().size();
    for (int i = 0; i < numParticles; i++) {
        ASSERT_EQUAL_VEC(s1.getPositions()[i], s2.getPositions()[i], TOL);
        ASSERT_EQUAL_VEC(s1.getVelocities()[i], s2.getVelocities()[i], TOL);
        Vec3 a1, b1, c1, a2, b2, c2;
        s1.getPeriodicBoxVectors(a1, b1, c1);
        s2.getPeriodicBoxVectors(a2, b2, c2);
        ASSERT_EQUAL_VEC(a1, a2, TOL);
        ASSERT_EQUAL_VEC(b1, b2, TOL);
        ASSERT_EQUAL_VEC(c1, c2, TOL);
        for (map<string, double>::const_iterator iter = s1.getParameters().begin(); iter != s1.getParameters().end(); ++iter)
68
            ASSERT_EQUAL(iter->second, (*s2.getParameters().find(iter->first)).second);
Peter Eastman's avatar
Peter Eastman committed
69
70
71
72
    }
}

void testCheckpoint() {
Peter Eastman's avatar
Peter Eastman committed
73
74
    const int numParticles = 100;
    const double boxSize = 5.0;
Peter Eastman's avatar
Peter Eastman committed
75
76
77
78
79
80
81
82
83
84
85
86
    const double temperature = 200.0;
    System system;
    system.addForce(new AndersenThermostat(0.0, 100.0));
    NonbondedForce* nonbonded = new NonbondedForce();
    system.addForce(nonbonded);
    nonbonded->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
    vector<Vec3> positions(numParticles);
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    for (int i = 0; i < numParticles; i++) {
        system.addParticle(1.0);
        nonbonded->addParticle(i%2 == 0 ? 0.1 : -0.1, 0.2, 0.1);
Peter Eastman's avatar
Peter Eastman committed
87
88
89
90
91
92
93
94
95
96
        bool clash;
        do {
            clash = false;
            positions[i] = Vec3(boxSize*genrand_real2(sfmt), boxSize*genrand_real2(sfmt), boxSize*genrand_real2(sfmt));
            for (int j = 0; j < i; j++) {
                Vec3 delta = positions[i]-positions[j];
                if (sqrt(delta.dot(delta)) < 0.1)
                    clash = true;
            }
        } while (clash);
Peter Eastman's avatar
Peter Eastman committed
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
    }
    VerletIntegrator integrator(0.001);
    Context context(system, integrator, platform);
    context.setPositions(positions);
    context.setPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
    context.setParameter(AndersenThermostat::Temperature(), temperature);
    
    // Run for a little while.
    
    integrator.step(100);
    
    // Record the current state and make a checkpoint.
    
    State s1 = context.getState(State::Positions | State::Velocities | State::Parameters);
    stringstream stream1(ios_base::out | ios_base::in | ios_base::binary);
    context.createCheckpoint(stream1);
    
    // Continue the simulation for a few more steps and record the state again.
    
    integrator.step(10);
    State s2 = context.getState(State::Positions | State::Velocities | State::Parameters);
    
    // Restore from the checkpoint and see if everything gets restored correctly.
    
    context.setPeriodicBoxVectors(Vec3(2*boxSize, 0, 0), Vec3(0, 2*boxSize, 0), Vec3(0, 0, 2*boxSize));
    context.setParameter(AndersenThermostat::Temperature(), temperature+10);
    context.loadCheckpoint(stream1);
    State s3 = context.getState(State::Positions | State::Velocities | State::Parameters);
    compareStates(s1, s3);
    
    // Now simulate from there and see if the trajectory is identical.
    
    integrator.step(10);
    State s4 = context.getState(State::Positions | State::Velocities | State::Parameters);
    compareStates(s2, s4);
Peter Eastman's avatar
Peter Eastman committed
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
159
    
    // Create a new Context that uses multiple devices.

    string deviceIndex = platform.getPropertyValue(context, OpenCLPlatform::OpenCLDeviceIndex());
    map<string, string> props;
    props[OpenCLPlatform::OpenCLDeviceIndex()] = deviceIndex+","+deviceIndex;
    VerletIntegrator integrator2(0.001);
    Context context2(system, integrator2, platform, props);
    context2.setPositions(positions);
    context2.setPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
    context2.setParameter(AndersenThermostat::Temperature(), temperature);
    
    // Now repeat all of the above tests with it.

    integrator2.step(100);
    State s5 = context2.getState(State::Positions | State::Velocities | State::Parameters);
    stringstream stream2(ios_base::out | ios_base::in | ios_base::binary);
    context2.createCheckpoint(stream2);
    integrator2.step(10);
    State s6 = context2.getState(State::Positions | State::Velocities | State::Parameters);
    context2.setPeriodicBoxVectors(Vec3(2*boxSize, 0, 0), Vec3(0, 2*boxSize, 0), Vec3(0, 0, 2*boxSize));
    context2.setParameter(AndersenThermostat::Temperature(), temperature+10);
    context2.loadCheckpoint(stream2);
    State s7 = context2.getState(State::Positions | State::Velocities | State::Parameters);
    compareStates(s5, s7);
    integrator2.step(10);
    State s8 = context2.getState(State::Positions | State::Velocities | State::Parameters);
    compareStates(s6, s8);
Peter Eastman's avatar
Peter Eastman committed
160
161
}

Peter Eastman's avatar
Peter Eastman committed
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
void testSetState() {
    const int numParticles = 10;
    const double boxSize = 3.0;
    const double temperature = 200.0;
    System system;
    system.addForce(new AndersenThermostat(0.0, 100.0));
    NonbondedForce* nonbonded = new NonbondedForce();
    system.addForce(nonbonded);
    nonbonded->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
    vector<Vec3> positions(numParticles);
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    for (int i = 0; i < numParticles; i++) {
        system.addParticle(1.0);
        nonbonded->addParticle(i%2 == 0 ? 0.1 : -0.1, 0.2, 0.1);
        positions[i] = Vec3(boxSize*genrand_real2(sfmt), boxSize*genrand_real2(sfmt), boxSize*genrand_real2(sfmt));
    }
    VerletIntegrator integrator(0.001);
    Context context(system, integrator, platform);
    context.setPositions(positions);
    context.setPeriodicBoxVectors(Vec3(boxSize, 0, 0), Vec3(0, boxSize, 0), Vec3(0, 0, boxSize));
    context.setParameter(AndersenThermostat::Temperature(), temperature);
    
    // Run for a little while.
    
    integrator.step(100);
    
    // Record the current state.
    
    State s1 = context.getState(State::Positions | State::Velocities | State::Parameters);
    
    // Continue the simulation for a few more steps and record a partial state.
    
    integrator.step(10);
    State s2 = context.getState(State::Positions);
    
    // Restore the original state and see if everything gets restored correctly.
    
    context.setPeriodicBoxVectors(Vec3(2*boxSize, 0, 0), Vec3(0, 2*boxSize, 0), Vec3(0, 0, 2*boxSize));
    context.setParameter(AndersenThermostat::Temperature(), temperature+10);
    context.setState(s1);
    State s3 = context.getState(State::Positions | State::Velocities | State::Parameters);
    compareStates(s1, s3);
    
    // Set the partial state and see if the correct things were set.
    
    context.setState(s2);
    State s4 = context.getState(State::Positions | State::Velocities | State::Parameters);
    for (int i = 0; i < numParticles; i++) {
        ASSERT_EQUAL_VEC(s2.getPositions()[i], s4.getPositions()[i], TOL);
        ASSERT_EQUAL_VEC(s1.getVelocities()[i], s4.getVelocities()[i], TOL);
    }
}

216
int main(int argc, char* argv[]) {
Peter Eastman's avatar
Peter Eastman committed
217
    try {
218
219
        if (argc > 1)
            platform.setPropertyDefaultValue("OpenCLPrecision", string(argv[1]));
Peter Eastman's avatar
Peter Eastman committed
220
        testCheckpoint();
Peter Eastman's avatar
Peter Eastman committed
221
        testSetState();
Peter Eastman's avatar
Peter Eastman committed
222
223
224
225
226
227
228
229
    }
    catch(const exception& e) {
        cout << "exception: " << e.what() << endl;
        return 1;
    }
    cout << "Done" << endl;
    return 0;
}