/* -------------------------------------------------------------------------- * * 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 reference implementation of CustomBondForce. */ #include "../../../tests/AssertionUtilities.h" #include "openmm/Context.h" #include "ReferencePlatform.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() { ReferencePlatform 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); } int main() { try { testBonds(); } catch(const exception& e) { cout << "exception: " << e.what() << endl; return 1; } cout << "Done" << endl; return 0; }