void testSerialization() {
    // Create a Force.

    GBVIForce force;
    force.setNonbondedMethod(GBVIForce::CutoffPeriodic);
    force.setBornRadiusScalingMethod(GBVIForce::QuinticSpline);
    force.setQuinticLowerLimitFactor(0.123);
    force.setQuinticUpperBornRadiusLimit(5.123);
    force.setCutoffDistance(2.0);
    force.setSoluteDielectric(5.1);
    force.setSolventDielectric(50.0);
    force.addParticle(1, 0.1, 0.01);
    force.addParticle(0.5, 0.2, 0.02);
    force.addParticle(-0.5, 0.3, 0.03);
    force.addBond(0, 1, 2.0);
    force.addBond(3, 5, 1.2);

    // Serialize and then deserialize it.

    stringstream buffer;
    XmlSerializer::serialize<GBVIForce>(&force, "Force", buffer);
    GBVIForce* copy = XmlSerializer::deserialize<GBVIForce>(buffer);

    // Compare the two forces to see if they are identical.

    GBVIForce& force2 = *copy;
    ASSERT_EQUAL(force.getNonbondedMethod(), force2.getNonbondedMethod());
    ASSERT_EQUAL(force.getCutoffDistance(), force2.getCutoffDistance());
    ASSERT_EQUAL(force.getSoluteDielectric(), force2.getSoluteDielectric());
    ASSERT_EQUAL(force.getSolventDielectric(), force2.getSolventDielectric());
    ASSERT_EQUAL(force.getNumParticles(), force2.getNumParticles());
    ASSERT_EQUAL(force.getQuinticUpperBornRadiusLimit(), force2.getQuinticUpperBornRadiusLimit());
    ASSERT_EQUAL(force.getQuinticLowerLimitFactor(), force2.getQuinticLowerLimitFactor());
    ASSERT_EQUAL(force.getBornRadiusScalingMethod(), force2.getBornRadiusScalingMethod());
    for (int i = 0; i < force.getNumParticles(); i++) {
        double charge1, radius1, scale1;
        double charge2, radius2, scale2;
        force.getParticleParameters(i, charge1, radius1, scale1);
        force2.getParticleParameters(i, charge2, radius2, scale2);
        ASSERT_EQUAL(charge1, charge2);
        ASSERT_EQUAL(radius1, radius2);
        ASSERT_EQUAL(scale1, scale2);
    }
    ASSERT_EQUAL(force.getNumBonds(), force2.getNumBonds());
    for (int i = 0; i < force.getNumBonds(); i++) {
        int a1, a2, b1, b2;
        double da, db;
        force.getBondParameters(i, a1, a2, da);
        force2.getBondParameters(i, b1, b2, db);
        ASSERT_EQUAL(a1, b1);
        ASSERT_EQUAL(a2, b2);
        ASSERT_EQUAL(da, db);
    }
}
Exemplo n.º 2
0
void testSingleParticle() {
    ReferencePlatform platform;
    System system;
    system.addParticle(2.0);
    LangevinIntegrator integrator(0, 0.1, 0.01);

    GBVIForce* forceField = new GBVIForce();

    double charge         = -1.0;
    double radius         =  0.15;
    double gamma          =  1.0;
    forceField->addParticle(charge, radius, gamma);
    system.addForce(forceField);
    ASSERT(!forceField->usesPeriodicBoundaryConditions());
    ASSERT(!system.usesPeriodicBoundaryConditions());

    Context context(system, integrator, platform);
    vector<Vec3> positions(1);
    positions[0] = Vec3(0, 0, 0);
    context.setPositions(positions);
    State state = context.getState(State::Energy);

    double bornRadius     = radius; 
    double eps0           = EPSILON0;
    double tau            = (1.0/forceField->getSoluteDielectric()-1.0/forceField->getSolventDielectric());

    double bornEnergy     = (-charge*charge/(8*PI_M*eps0))*tau/bornRadius;
    double nonpolarEnergy = -gamma*tau*std::pow(radius/bornRadius, 3.0);

    double expectedE      = (bornEnergy+nonpolarEnergy); 
    double obtainedE      = state.getPotentialEnergy(); 
    double diff           = fabs((obtainedE - expectedE)/expectedE);

    ASSERT_EQUAL_TOL((bornEnergy+nonpolarEnergy), state.getPotentialEnergy(), 0.01);
}
Exemplo n.º 3
0
void testEnergyEthane(int applyBornRadiiScaling) {

    ReferencePlatform platform;
    const int numParticles = 8;
    System system;
    LangevinIntegrator integrator(0, 0.1, 0.01);

    // harmonic bond

    double C_HBondDistance   = 0.1097;
    double C_CBondDistance   = 0.1504;
    HarmonicBondForce* bonds = new HarmonicBondForce();
    bonds->addBond(0, 1, C_HBondDistance, 0.0);
    bonds->addBond(2, 1, C_HBondDistance, 0.0);
    bonds->addBond(3, 1, C_HBondDistance, 0.0);

    bonds->addBond(1, 4, C_CBondDistance, 0.0);

    bonds->addBond(5, 4, C_HBondDistance, 0.0);
    bonds->addBond(6, 4, C_HBondDistance, 0.0);
    bonds->addBond(7, 4, C_HBondDistance, 0.0);

    system.addForce(bonds);

    double C_radius, C_gamma, C_charge, H_radius, H_gamma, H_charge;

    int AM1_BCC = 1;
    H_charge    = -0.053;
    C_charge    = -3.0*H_charge;
    if (AM1_BCC) {
       C_radius =  0.180;
       C_gamma  = -0.2863;
       H_radius =  0.125;
       H_gamma  =  0.2437;
    }
    else {
       C_radius =  0.215;
       C_gamma  = -1.1087;
       H_radius =  0.150;
       H_gamma  =  0.1237;
    }

    NonbondedForce* nonbonded = new NonbondedForce();
    nonbonded->setNonbondedMethod(NonbondedForce::NoCutoff);

    GBVIForce* forceField = new GBVIForce();
    if (applyBornRadiiScaling) {
        forceField->setBornRadiusScalingMethod(GBVIForce::QuinticSpline);
    }
    else {
        forceField->setBornRadiusScalingMethod(GBVIForce::NoScaling);
    }
    for (int i = 0; i < numParticles; i++) {
       system.addParticle(1.0);
       forceField->addParticle(H_charge, H_radius, H_gamma);
       nonbonded->addParticle( H_charge, H_radius, 0.0);
    }
 
    forceField->setParticleParameters(1, C_charge, C_radius, C_gamma);
    forceField->setParticleParameters(4, C_charge, C_radius, C_gamma);
 
    nonbonded->setParticleParameters( 1, C_charge, C_radius, 0.0);
    nonbonded->setParticleParameters( 4, C_charge, C_radius, 0.0);
 
    forceField->addBond(0, 1, C_HBondDistance);
    forceField->addBond(2, 1, C_HBondDistance);
    forceField->addBond(3, 1, C_HBondDistance);
    forceField->addBond(1, 4, C_CBondDistance);
    forceField->addBond(5, 4, C_HBondDistance);
    forceField->addBond(6, 4, C_HBondDistance);
    forceField->addBond(7, 4, C_HBondDistance);
    
    std::vector<pair<int, int> > bondExceptions;
    std::vector<double> bondDistances;
    
    bondExceptions.push_back(pair<int, int>(0, 1)); 
    bondDistances.push_back(C_HBondDistance);
    
    bondExceptions.push_back(pair<int, int>(2, 1)); 
    bondDistances.push_back(C_HBondDistance);
    
    bondExceptions.push_back(pair<int, int>(3, 1)); 
    bondDistances.push_back(C_HBondDistance);
    
    bondExceptions.push_back(pair<int, int>(1, 4)); 
    bondDistances.push_back(C_CBondDistance);
    
    bondExceptions.push_back(pair<int, int>(5, 4)); 
    bondDistances.push_back(C_HBondDistance);
    
    bondExceptions.push_back(pair<int, int>(6, 4)); 
    bondDistances.push_back(C_HBondDistance);
 
    bondExceptions.push_back(pair<int, int>(7, 4));
    bondDistances.push_back(C_HBondDistance);
 
    nonbonded->createExceptionsFromBonds(bondExceptions, 0.0, 0.0);
 
    system.addForce(forceField);
    system.addForce(nonbonded);

    Context context(system, integrator, platform);
    
    vector<Vec3> positions(numParticles);
    positions[0] = Vec3(0.5480,    1.7661,    0.0000);
    positions[1] = Vec3(0.7286,    0.8978,    0.6468);
    positions[2] = Vec3(0.4974,    0.0000,    0.0588);
    positions[3] = Vec3(0.0000,    0.9459,    1.4666);
    positions[4] = Vec3(2.1421,    0.8746,    1.1615);
    positions[5] = Vec3(2.3239,    0.0050,    1.8065);
    positions[6] = Vec3(2.8705,    0.8295,    0.3416);
    positions[7] = Vec3(2.3722,    1.7711,    1.7518);
    context.setPositions(positions);

    State state = context.getState(State::Forces | State::Energy);
    
    // Take a small step in the direction of the energy gradient.
    
    double norm        = 0.0;
    double forceSum[3] = { 0.0, 0.0, 0.0 };
    for (int i = 0; i < numParticles; ++i) {
        Vec3 f  = state.getForces()[i];
        norm        += f[0]*f[0] + f[1]*f[1] + f[2]*f[2];
        forceSum[0] += f[0];
        forceSum[1] += f[1];
        forceSum[2] += f[2];
    }
    norm               = std::sqrt(norm);

    const double delta = 1e-4;
    double step = delta/norm;
    for (int i = 0; i < numParticles; ++i) {
        Vec3 p = positions[i];
        Vec3 f = state.getForces()[i];
        positions[i] = Vec3(p[0]-f[0]*step, p[1]-f[1]*step, p[2]-f[2]*step);
    }
    context.setPositions(positions);
    
    State state2 = context.getState(State::Energy);

    // See whether the potential energy changed by the expected amount.
    
    ASSERT_EQUAL_TOL(norm, (state2.getPotentialEnergy()-state.getPotentialEnergy())/delta, 0.01)
}