void testIdealGas() { const int numParticles = 64; const int frequency = 10; const int steps = 1000; const double pressure = 1.5; const double pressureInMD = pressure*(AVOGADRO*1e-25); // pressure in kJ/mol/nm^3 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)); } MonteCarloAnisotropicBarostat* barostat = new MonteCarloAnisotropicBarostat(Vec3(pressure, pressure, pressure), temp[0], true, true, true, frequency); 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]; context.getState(0).getPeriodicBoxVectors(box[0], box[1], box[2]); volume += box[0][0]*box[1][1]*box[2][2]; integrator.step(frequency); } volume /= steps; double expected = (numParticles+1)*BOLTZ*temp[i]/pressureInMD; ASSERT_USUALLY_EQUAL_TOL(expected, volume, 3/std::sqrt((double) steps)); } }
/** * Run a constant pressure simulation on an anisotropic Einstein crystal * using isotropic and anisotropic barostats. There are a total of 15 simulations: * * 1) 3 pressures: 9.0, 10.0, 11.0 bar, for each of the following groups: * 2) 3 groups of simulations that scale just one axis: x, y, z * 3) 1 group of simulations that scales all three axes in the anisotropic barostat * 4) 1 group of simulations that scales all three axes in the isotropic barostat * * Results that we will check: * * a) In each group of simulations, the volume should decrease with increasing pressure * b) In the three simulation groups that scale just one axis, the compressibility (i.e. incremental volume change * with increasing pressure) should go like kx > ky > kz (because the spring constant is largest in the z-direction) * c) The anisotropic barostat should produce the same result as the isotropic barostat when all three axes are scaled */ void testEinsteinCrystal() { const int numParticles = 64; const int frequency = 2; const int equil = 10000; const int steps = 5000; const double pressure = 10.0; const double pressureInMD = pressure*(AVOGADRO*1e-25); // pressure in kJ/mol/nm^3 const double temp = 300.0; // Only test one temperature since we're looking at three pressures. const double pres3[] = {2.0, 8.0, 15.0}; const double initialVolume = numParticles*BOLTZ*temp/pressureInMD; const double initialLength = std::pow(initialVolume, 1.0/3.0); vector<double> initialPositions(3); vector<double> results; // Run four groups of anisotropic simulations; scaling just x, y, z, then all three. for (int a = 0; a < 4; a++) { // Test barostat for three different pressures. for (int p = 0; p < 3; p++) { // Create a system of noninteracting particles attached by harmonic springs to their initial positions. System system; system.setDefaultPeriodicBoxVectors(Vec3(initialLength, 0, 0), Vec3(0, initialLength, 0), Vec3(0, 0, initialLength)); vector<Vec3> positions(numParticles); OpenMM_SFMT::SFMT sfmt; init_gen_rand(0, sfmt); // Anisotropic force constants. CustomExternalForce* force = new CustomExternalForce("0.005*(x-x0)^2 + 0.01*(y-y0)^2 + 0.02*(z-z0)^2"); force->addPerParticleParameter("x0"); force->addPerParticleParameter("y0"); force->addPerParticleParameter("z0"); NonbondedForce* nb = new NonbondedForce(); nb->setNonbondedMethod(NonbondedForce::CutoffPeriodic); for (int i = 0; i < numParticles; ++i) { system.addParticle(1.0); positions[i] = Vec3(((i/16)%4+0.5)*initialLength/4, ((i/4)%4+0.5)*initialLength/4, (i%4+0.5)*initialLength/4); initialPositions[0] = positions[i][0]; initialPositions[1] = positions[i][1]; initialPositions[2] = positions[i][2]; force->addParticle(i, initialPositions); nb->addParticle(0, initialLength/6, 0.1); } system.addForce(force); system.addForce(nb); // Create the barostat. MonteCarloAnisotropicBarostat* barostat = new MonteCarloAnisotropicBarostat(Vec3(pres3[p], pres3[p], pres3[p]), temp, (a==0||a==3), (a==1||a==3), (a==2||a==3), frequency); system.addForce(barostat); barostat->setTemperature(temp); LangevinIntegrator integrator(temp, 0.1, 0.01); Context context(system, integrator, platform); context.setPositions(positions); // Let it equilibrate. integrator.step(equil); // 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]; context.getState(0).getPeriodicBoxVectors(box[0], box[1], box[2]); volume += box[0][0]*box[1][1]*box[2][2]; integrator.step(frequency); } volume /= steps; results.push_back(volume); } } for (int p = 0; p < 3; p++) { // Create a system of noninteracting particles attached by harmonic springs to their initial positions. System system; system.setDefaultPeriodicBoxVectors(Vec3(initialLength, 0, 0), Vec3(0, initialLength, 0), Vec3(0, 0, initialLength)); vector<Vec3> positions(numParticles); OpenMM_SFMT::SFMT sfmt; init_gen_rand(0, sfmt); // Anisotropic force constants. CustomExternalForce* force = new CustomExternalForce("0.005*(x-x0)^2 + 0.01*(y-y0)^2 + 0.02*(z-z0)^2"); force->addPerParticleParameter("x0"); force->addPerParticleParameter("y0"); force->addPerParticleParameter("z0"); NonbondedForce* nb = new NonbondedForce(); nb->setNonbondedMethod(NonbondedForce::CutoffPeriodic); for (int i = 0; i < numParticles; ++i) { system.addParticle(1.0); positions[i] = Vec3(((i/16)%4+0.5)*initialLength/4, ((i/4)%4+0.5)*initialLength/4, (i%4+0.5)*initialLength/4); initialPositions[0] = positions[i][0]; initialPositions[1] = positions[i][1]; initialPositions[2] = positions[i][2]; force->addParticle(i, initialPositions); nb->addParticle(0, initialLength/6, 0.1); } system.addForce(force); system.addForce(nb); // Create the barostat. MonteCarloBarostat* barostat = new MonteCarloBarostat(pres3[p], temp, frequency); system.addForce(barostat); barostat->setTemperature(temp); LangevinIntegrator integrator(temp, 0.1, 0.001); Context context(system, integrator, platform); context.setPositions(positions); // Let it equilibrate. integrator.step(equil); // 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]; context.getState(0).getPeriodicBoxVectors(box[0], box[1], box[2]); volume += box[0][0]*box[1][1]*box[2][2]; integrator.step(frequency); } volume /= steps; results.push_back(volume); } // Check to see if volumes decrease with increasing pressure. ASSERT_USUALLY_TRUE(results[0] > results[1]); ASSERT_USUALLY_TRUE(results[1] > results[2]); ASSERT_USUALLY_TRUE(results[3] > results[4]); ASSERT_USUALLY_TRUE(results[4] > results[5]); ASSERT_USUALLY_TRUE(results[6] > results[7]); ASSERT_USUALLY_TRUE(results[7] > results[8]); // Check to see if incremental volume changes with increasing pressure go like kx > ky > kz. ASSERT_USUALLY_TRUE((results[0] - results[1]) > (results[3] - results[4])); ASSERT_USUALLY_TRUE((results[1] - results[2]) > (results[4] - results[5])); ASSERT_USUALLY_TRUE((results[3] - results[4]) > (results[6] - results[7])); ASSERT_USUALLY_TRUE((results[4] - results[5]) > (results[7] - results[8])); // Check to see if the volumes are equal for isotropic and anisotropic (all axis). ASSERT_USUALLY_EQUAL_TOL(results[9], results[12], 3/std::sqrt((double) steps)); ASSERT_USUALLY_EQUAL_TOL(results[10], results[13], 3/std::sqrt((double) steps)); ASSERT_USUALLY_EQUAL_TOL(results[11], results[14], 3/std::sqrt((double) steps)); }