/* -------------------------------------------------------------------------- * * 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-2009 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 OpenCL implementation of CustomBondForce. */ #include "openmm/internal/AssertionUtilities.h" #include "openmm/Context.h" #include "OpenCLPlatform.h" #include "openmm/CustomBondForce.h" #include "openmm/System.h" #include "openmm/VerletIntegrator.h" #include "../src/SimTKUtilities/SimTKOpenMMRealType.h" #include #include using namespace OpenMM; using namespace std; const double TOL = 1e-5; void testBonds() { OpenCLPlatform platform; System system; system.addParticle(1.0); system.addParticle(1.0); system.addParticle(1.0); VerletIntegrator integrator(0.01); CustomBondForce* forceField = new CustomBondForce("scale*k*(r-r0)^2"); forceField->addPerBondParameter("r0"); forceField->addPerBondParameter("k"); forceField->addGlobalParameter("scale", 0.5); vector parameters(2); parameters[0] = 1.5; parameters[1] = 0.8; forceField->addBond(0, 1, parameters); parameters[0] = 1.2; parameters[1] = 0.7; forceField->addBond(1, 2, parameters); system.addForce(forceField); Context context(system, integrator, platform); vector positions(3); positions[0] = Vec3(0, 2, 0); positions[1] = Vec3(0, 0, 0); positions[2] = Vec3(1, 0, 0); context.setPositions(positions); State state = context.getState(State::Forces | State::Energy); { const vector& forces = state.getForces(); ASSERT_EQUAL_VEC(Vec3(0, -0.8*0.5, 0), forces[0], TOL); ASSERT_EQUAL_VEC(Vec3(0.7*0.2, 0, 0), forces[2], TOL); ASSERT_EQUAL_VEC(Vec3(-forces[0][0]-forces[2][0], -forces[0][1]-forces[2][1], -forces[0][2]-forces[2][2]), forces[1], TOL); ASSERT_EQUAL_TOL(0.5*0.8*0.5*0.5 + 0.5*0.7*0.2*0.2, state.getPotentialEnergy(), TOL); } // Try changing the bond parameters and make sure it's still correct. parameters[0] = 1.6; parameters[1] = 0.9; forceField->setBondParameters(0, 0, 1, parameters); parameters[0] = 1.3; parameters[1] = 0.8; forceField->setBondParameters(1, 1, 2, parameters); forceField->updateParametersInContext(context); state = context.getState(State::Forces | State::Energy); { const vector& forces = state.getForces(); ASSERT_EQUAL_VEC(Vec3(0, -0.9*0.4, 0), forces[0], TOL); ASSERT_EQUAL_VEC(Vec3(0.8*0.3, 0, 0), forces[2], TOL); ASSERT_EQUAL_VEC(Vec3(-forces[0][0]-forces[2][0], -forces[0][1]-forces[2][1], -forces[0][2]-forces[2][2]), forces[1], TOL); ASSERT_EQUAL_TOL(0.5*0.9*0.4*0.4 + 0.5*0.8*0.3*0.3, state.getPotentialEnergy(), TOL); } } void testManyParameters() { OpenCLPlatform platform; System system; system.addParticle(1.0); system.addParticle(1.0); VerletIntegrator integrator(0.01); CustomBondForce* forceField = new CustomBondForce("(a+b+c+d+e+f+g+h+i)*r"); forceField->addPerBondParameter("a"); forceField->addPerBondParameter("b"); forceField->addPerBondParameter("c"); forceField->addPerBondParameter("d"); forceField->addPerBondParameter("e"); forceField->addPerBondParameter("f"); forceField->addPerBondParameter("g"); forceField->addPerBondParameter("h"); forceField->addPerBondParameter("i"); vector parameters(forceField->getNumPerBondParameters()); for (int i = 0; i < parameters.size(); i++) parameters[i] = i; forceField->addBond(0, 1, parameters); system.addForce(forceField); Context context(system, integrator, platform); vector positions(2); positions[0] = Vec3(0, 0, 0); positions[1] = Vec3(0, 2.5, 0); context.setPositions(positions); State state = context.getState(State::Forces | State::Energy); const vector& forces = state.getForces(); double f = 1+2+3+4+5+6+7+8; ASSERT_EQUAL_VEC(Vec3(0, f, 0), forces[0], TOL); ASSERT_EQUAL_VEC(Vec3(0, -f, 0), forces[1], TOL); ASSERT_EQUAL_TOL(f*2.5, state.getPotentialEnergy(), TOL); } void testParallelComputation() { OpenCLPlatform platform; System system; const int numParticles = 200; for (int i = 0; i < numParticles; i++) system.addParticle(1.0); CustomBondForce* force = new CustomBondForce(("(r-1.1)^2")); vector params; for (int i = 1; i < numParticles; i++) force->addBond(i-1, i, params); system.addForce(force); vector positions(numParticles); for (int i = 0; i < numParticles; i++) positions[i] = Vec3(i, 0, 0); VerletIntegrator integrator1(0.01); Context context1(system, integrator1, platform); context1.setPositions(positions); State state1 = context1.getState(State::Forces | State::Energy); VerletIntegrator integrator2(0.01); string deviceIndex = platform.getPropertyValue(context1, OpenCLPlatform::OpenCLDeviceIndex()); map props; props[OpenCLPlatform::OpenCLDeviceIndex()] = deviceIndex+","+deviceIndex; Context context2(system, integrator2, platform, props); context2.setPositions(positions); State state2 = context2.getState(State::Forces | State::Energy); ASSERT_EQUAL_TOL(state1.getPotentialEnergy(), state2.getPotentialEnergy(), 1e-5); for (int i = 0; i < numParticles; i++) ASSERT_EQUAL_VEC(state1.getForces()[i], state2.getForces()[i], 1e-5); } int main() { try { testBonds(); testManyParameters(); testParallelComputation(); } catch(const exception& e) { cout << "exception: " << e.what() << endl; return 1; } cout << "Done" << endl; return 0; }