HelloArgon.cpp 2.59 KB
Newer Older
1
2
3
#include "OpenMM.h"
#include <cstdio>

4
// forward declaration of subroutine defined later in this file.
5
void writePdb(const OpenMM::OpenMMContext& context);
6
7
8
9

void simulateArgon() 
{
    // Load any shared libraries containing GPU implementations
10
11
    OpenMM::Platform::loadPluginsFromDirectory(
        OpenMM::Platform::getDefaultPluginsDirectory());
12

13
14
    OpenMM::System system;
    OpenMM::NonbondedForce* nonbond = new OpenMM::NonbondedForce(); 
15
16
    system.addForce(nonbond);

17
    // Create three atoms
18
    std::vector<OpenMM::Vec3> initialPositions(3);
19
20
    for (int a = 0; a < 3; ++a) 
    {
21
        system.addParticle(39.95); // mass
22
23
24
25
26

        // charge, sigma, well depth
        nonbond->addParticle(0.0, 0.3350, 0.001603); 

        initialPositions[a] = OpenMM::Vec3(0.5*a,0,0); // location
27
28
    }

29
    OpenMM::VerletIntegrator integrator(0.020); // step size in picoseconds
30

31
    // Let OpenMM Context choose best platform.
32
    OpenMM::OpenMMContext context(system, integrator);
33
34
    printf( "REMARK  Using OpenMM platform %s\n", 
        context.getPlatform().getName().c_str() );
35

36
37
    // Set the starting positions of the atoms. Velocities will be zero.
    context.setPositions(initialPositions);
38

39
    // Simulate
40
    while(context.getTime() < 500.0) { // picoseconds
41
42
43
44
45
46
47
48
49
50
51
52
        writePdb(context); // output coordinates
        // Run 100 steps at a time, for efficient use of OpenMM
        integrator.step(100);
    }
    writePdb(context); // output final coordinates
}

int main() 
{
    try {
        simulateArgon();
        return 0; // success!
53
    }
54
55
56
57
58
59
60
61
    // Catch and report usage and runtime errors detected by OpenMM and fail.
    catch(const std::exception& e) {
        printf("EXCEPTION: %s\n", e.what());
        return 1; // failure!
    }
}

// writePdb() subroutine for quick-and-dirty trajectory output.
62
void writePdb(const OpenMM::OpenMMContext& context) 
63
64
{
    // Request atomic positions from OpenMM
65
66
    const OpenMM::State state = context.getState(OpenMM::State::Positions);
    const std::vector<OpenMM::Vec3>& pos = state.getPositions();
67

68
69
70
71
72
73
74
75
76
77
78
79
80
81
    // write out in PDB format

    // Use PDB MODEL cards to number trajectory frames
    static int modelFrameNumber = 0; 
    modelFrameNumber++;
    printf("MODEL     %d\n", modelFrameNumber); // start of frame
    for (int a = 0; a < context.getSystem().getNumParticles(); ++a)
    {
        printf("ATOM  %5d AR    AR     1    ", a+1); // atom number
        printf("%8.3f%8.3f%8.3f  1.00  0.00          AR\n",
            // notice "*10" converts nanometers to Angstroms
            pos[a][0]*10, pos[a][1]*10, pos[a][2]*10);
    }
    printf("ENDMDL\n"); // end of frame
82
}