/* -------------------------------------------------------------------------- * * 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-2015 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. * * -------------------------------------------------------------------------- */ #include "CudaTests.h" #include "TestNonbondedForce.h" void testParallelComputation(NonbondedForce::NonbondedMethod method) { System system; const int numParticles = 200; for (int i = 0; i < numParticles; i++) system.addParticle(1.0); NonbondedForce* force = new NonbondedForce(); for (int i = 0; i < numParticles; i++) force->addParticle(i%2-0.5, 0.5, 1.0); force->setNonbondedMethod(method); system.addForce(force); system.setDefaultPeriodicBoxVectors(Vec3(5,0,0), Vec3(0,5,0), Vec3(0,0,5)); OpenMM_SFMT::SFMT sfmt; init_gen_rand(0, sfmt); vector positions(numParticles); for (int i = 0; i < numParticles; i++) positions[i] = Vec3(5*genrand_real2(sfmt), 5*genrand_real2(sfmt), 5*genrand_real2(sfmt)); for (int i = 0; i < numParticles; ++i) for (int j = 0; j < i; ++j) { Vec3 delta = positions[i]-positions[j]; if (delta.dot(delta) < 0.1) force->addException(i, j, 0, 1, 0); } // Create two contexts, one with a single device and one with two devices. 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, CudaPlatform::CudaDeviceIndex()); map props; props[CudaPlatform::CudaDeviceIndex()] = deviceIndex+","+deviceIndex; Context context2(system, integrator2, platform, props); context2.setPositions(positions); State state2 = context2.getState(State::Forces | State::Energy); // See if they agree. 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); // Modify some particle parameters and see if they still agree. for (int i = 0; i < numParticles; i += 5) { double charge, sigma, epsilon; force->getParticleParameters(i, charge, sigma, epsilon); force->setParticleParameters(i, 0.9*charge, sigma, epsilon); } force->updateParametersInContext(context1); force->updateParametersInContext(context2); state1 = context1.getState(State::Forces | State::Energy); 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); } void testReordering() { // Check that reordering of atoms doesn't alter their positions. const int numParticles = 200; System system; system.setDefaultPeriodicBoxVectors(Vec3(6, 0, 0), Vec3(2.1, 6, 0), Vec3(-1.5, -0.5, 6)); NonbondedForce *nonbonded = new NonbondedForce(); nonbonded->setNonbondedMethod(NonbondedForce::PME); system.addForce(nonbonded); vector positions; OpenMM_SFMT::SFMT sfmt; init_gen_rand(0, sfmt); for (int i = 0; i < numParticles; i++) { system.addParticle(1.0); nonbonded->addParticle(0.0, 0.0, 0.0); positions.push_back(Vec3(genrand_real2(sfmt)-0.5, genrand_real2(sfmt)-0.5, genrand_real2(sfmt)-0.5)*20); } VerletIntegrator integrator(0.001); Context context(system, integrator, platform); context.setPositions(positions); integrator.step(1); State state = context.getState(State::Positions | State::Velocities); for (int i = 0; i < numParticles; i++) { ASSERT_EQUAL_VEC(positions[i], state.getPositions()[i], 1e-6); } } void testDeterministicForces() { // Check that the CudaDeterministicForces property works correctly. const int numParticles = 1000; System system; system.setDefaultPeriodicBoxVectors(Vec3(6, 0, 0), Vec3(2.1, 6, 0), Vec3(-1.5, -0.5, 6)); NonbondedForce *nonbonded = new NonbondedForce(); nonbonded->setNonbondedMethod(NonbondedForce::PME); system.addForce(nonbonded); vector positions; 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 ? 1 : -1, 1, 0); positions.push_back(Vec3(genrand_real2(sfmt)-0.5, genrand_real2(sfmt)-0.5, genrand_real2(sfmt)-0.5)*6); } VerletIntegrator integrator(0.001); map properties; properties[CudaPlatform::CudaDeterministicForces()] = "true"; Context context(system, integrator, platform, properties); context.setPositions(positions); State state1 = context.getState(State::Forces); State state2 = context.getState(State::Forces); // All forces should be *exactly* equal. for (int i = 0; i < numParticles; i++) { ASSERT_EQUAL(state1.getForces()[i][0], state2.getForces()[i][0]); ASSERT_EQUAL(state1.getForces()[i][1], state2.getForces()[i][1]); ASSERT_EQUAL(state1.getForces()[i][2], state2.getForces()[i][2]); } } void runPlatformTests() { testParallelComputation(NonbondedForce::NoCutoff); testParallelComputation(NonbondedForce::Ewald); testParallelComputation(NonbondedForce::PME); testReordering(); testDeterministicForces(); }