"...git@developer.sourcefind.cn:2222/modelzoo/trellis.2.git" did not exist on "f05e915fd8beb0b629a2f7eb58a52304c47f6370"
TestReferenceMonteCarloBarostat.cpp 8.83 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-2010 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 MonteCarloBarostat.
 */

36
#include "openmm/internal/AssertionUtilities.h"
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
#include "openmm/MonteCarloBarostat.h"
#include "openmm/Context.h"
#include "ReferencePlatform.h"
#include "openmm/NonbondedForce.h"
#include "openmm/System.h"
#include "openmm/LangevinIntegrator.h"
#include "openmm/VerletIntegrator.h"
#include "sfmt/SFMT.h"
#include "../src/SimTKUtilities/SimTKOpenMMRealType.h"
#include <iostream>
#include <vector>

using namespace OpenMM;
using namespace std;

void testChangingBoxSize() {
    ReferencePlatform platform;
    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(4, 0, 0), Vec3(0, 5, 0), Vec3(0, 0, 6));
56
    system.addParticle(1.0);
57
58
59
60
61
    NonbondedForce* nb = new NonbondedForce();
    nb->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
    nb->setCutoffDistance(2.0);
    nb->addParticle(1, 0.5, 0.5);
    system.addForce(nb);
62
63
    LangevinIntegrator integrator(300.0, 1.0, 0.01);
    Context context(system, integrator, platform);
64
65
66
    vector<Vec3> positions;
    positions.push_back(Vec3());
    context.setPositions(positions);
67
    Vec3 x, y, z;
68
    context.getState(State::Forces).getPeriodicBoxVectors(x, y, z);
69
70
71
72
    ASSERT_EQUAL_VEC(Vec3(4, 0, 0), x, 0);
    ASSERT_EQUAL_VEC(Vec3(0, 5, 0), y, 0);
    ASSERT_EQUAL_VEC(Vec3(0, 0, 6), z, 0);
    context.setPeriodicBoxVectors(Vec3(7, 0, 0), Vec3(0, 8, 0), Vec3(0, 0, 9));
73
    context.getState(State::Forces).getPeriodicBoxVectors(x, y, z);
74
75
76
    ASSERT_EQUAL_VEC(Vec3(7, 0, 0), x, 0);
    ASSERT_EQUAL_VEC(Vec3(0, 8, 0), y, 0);
    ASSERT_EQUAL_VEC(Vec3(0, 0, 9), z, 0);
77
78
79
80
81
82
83
84
85
86
87
88
    
    // Shrinking the box too small should produce an exception.
    
    context.setPeriodicBoxVectors(Vec3(7, 0, 0), Vec3(0, 3.9, 0), Vec3(0, 0, 9));
    bool ok = true;
    try {
        context.getState(State::Forces).getPeriodicBoxVectors(x, y, z);
        ok = false;
    }
    catch (exception& ex) {
    }
    ASSERT(ok);
89
90
}

91
void testIdealGas(int aniso) {
92
93
94
95
    const int numParticles = 64;
    const int frequency = 10;
    const int steps = 1000;
    const double pressure = 1.5;
96
    const double pressureInMD = pressure*(AVOGADRO*1e-25); // pressure in kJ/mol/nm^3
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
    const double temp[] = {300.0, 600.0, 1000.0};
    const double initialVolume = numParticles*BOLTZ*temp[1]/pressureInMD;
    const double initialLength = std::pow(initialVolume, 1.0/3.0);

    // Create a gas of noninteracting particles.

    ReferencePlatform platform;
    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(initialLength, 0, 0), Vec3(0, 0.5*initialLength, 0), Vec3(0, 0, 2*initialLength));
    vector<Vec3> positions(numParticles);
    OpenMM_SFMT::SFMT sfmt;
    init_gen_rand(0, sfmt);
    for (int i = 0; i < numParticles; ++i) {
        system.addParticle(1.0);
        positions[i] = Vec3(initialLength*genrand_real2(sfmt), 0.5*initialLength*genrand_real2(sfmt), 2*initialLength*genrand_real2(sfmt));
    }
113
    MonteCarloBarostat* barostat = new MonteCarloBarostat(pressure, temp[0], frequency, aniso);
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
    system.addForce(barostat);

    // Test it for three different temperatures.

    for (int i = 0; i < 3; i++) {
        barostat->setTemperature(temp[i]);
        LangevinIntegrator integrator(temp[i], 0.1, 0.01);
        Context context(system, integrator, platform);
        context.setPositions(positions);

        // Let it equilibrate.

        integrator.step(10000);

        // Now run it for a while and see if the volume is correct.

        double volume = 0.0;
        for (int j = 0; j < steps; ++j) {
            Vec3 box[3];
133
            context.getState(0).getPeriodicBoxVectors(box[0], box[1], box[2]);
134
            volume += box[0][0]*box[1][1]*box[2][2];
135
136
137
138
	    if (!aniso) {
	      ASSERT_EQUAL_TOL(0.5*box[0][0], box[1][1], 1e-5);
	      ASSERT_EQUAL_TOL(2*box[0][0], box[2][2], 1e-5);
	    }
139
140
141
            integrator.step(frequency);
        }
        volume /= steps;
142
        double expected = (numParticles+1)*BOLTZ*temp[i]/pressureInMD;
143
        ASSERT_USUALLY_EQUAL_TOL(expected, volume, 3/std::sqrt((double) steps));
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
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
    }
}

void testRandomSeed() {
    const int numParticles = 8;
    const double temp = 100.0;
    const double pressure = 1.5;
    ReferencePlatform platform;
    System system;
    system.setDefaultPeriodicBoxVectors(Vec3(8, 0, 0), Vec3(0, 8, 0), Vec3(0, 0, 8));
    VerletIntegrator integrator(0.01);
    NonbondedForce* forceField = new NonbondedForce();
    forceField->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
    for (int i = 0; i < numParticles; ++i) {
        system.addParticle(2.0);
        forceField->addParticle((i%2 == 0 ? 1.0 : -1.0), 1.0, 5.0);
    }
    system.addForce(forceField);
    MonteCarloBarostat* barostat = new MonteCarloBarostat(pressure, temp, 1);
    system.addForce(barostat);
    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.

    barostat->setRandomNumberSeed(5);
    Context context(system, integrator, platform);
    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.

    barostat->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]);
        }
    }
}

int main() {
    try {
        testChangingBoxSize();
213
214
        testIdealGas(0);
        testIdealGas(1);
215
216
217
218
219
220
221
222
223
224
        testRandomSeed();
    }
    catch(const exception& e) {
        cout << "exception: " << e.what() << endl;
        return 1;
    }
    cout << "Done" << endl;
    return 0;
}