Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add compatibility with openforcefield 0.7.0 #121

Merged
merged 6 commits into from
Jun 24, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 3 additions & 0 deletions README.md
Original file line number Diff line number Diff line change
Expand Up @@ -269,6 +269,9 @@ See the corresponding directories for information on how to use the provided con

# Changelog

## 0.7.3 Bugfix release: Compatibility with openforcefield toolkit 0.7.0 and auto-detection of installed openforcefield force fields
* [(PR #121)](https://github.com/openmm/openmmforcefields/pull/121) Add compatibility with [`openforcefield 0.7.0`](https://github.com/openforcefield/openforcefield/releases/tag/0.7.0)

## 0.7.3 Bugfix release: Compatibility with openforcefield toolkit 0.7.0 and auto-detection of installed openforcefield force fields
* [(PR #119)](https://github.com/openmm/openmmforcefields/pull/119) Handle `None` partial charges in openforcefield `Molecule` objects (needed in `openforcefield` toolkit 0.7.0)
* [(PR #120)](https://github.com/openmm/openmmforcefields/pull/120) Auto-detect installed SMIRNOFF force fields
Expand Down
6 changes: 5 additions & 1 deletion openmmforcefields/generators/template_generators.py
Original file line number Diff line number Diff line change
Expand Up @@ -570,6 +570,10 @@ def generate_residue_template(self, molecule, residue_atoms=None):

# Compute net formal charge
net_charge = molecule.total_charge
from simtk import unit
if type(net_charge) != unit.Quantity:
# openforcefield toolkit < 0.7.0 did not return unit-bearing quantity
net_charge = float(net_charge) * unit.elementary_charge
_logger.debug(f'Total charge is {net_charge}')

# Compute partial charges if required
Expand Down Expand Up @@ -621,7 +625,7 @@ def generate_residue_template(self, molecule, residue_atoms=None):
residue_charge = 0.0 * unit.elementary_charge
total_charge = unit.sum(molecule.partial_charges)
sum_of_absolute_charge = unit.sum(abs(molecule.partial_charges))
charge_deficit = net_charge * unit.elementary_charge - total_charge
charge_deficit = net_charge - total_charge
if sum_of_absolute_charge / unit.elementary_charge > 0.0:
# Redistribute excess charge proportionally to absolute charge
molecule.partial_charges += charge_deficit * abs(molecule.partial_charges) / sum_of_absolute_charge
Expand Down
13 changes: 9 additions & 4 deletions openmmforcefields/tests/test_template_generators.py
Original file line number Diff line number Diff line change
Expand Up @@ -206,7 +206,7 @@ def test_charge(self):
from simtk import unit
molecule = self.molecules[0]
# Ensure partial charges are initially zero
assert np.all(molecule.partial_charges / unit.elementary_charge == 0)
assert (molecule.partial_charges is None) or np.all(molecule.partial_charges / unit.elementary_charge == 0)
# Add the molecule
generator.add_molecules(molecule)
# Create the System
Expand All @@ -231,9 +231,13 @@ def test_charge_from_molecules(self):
from simtk import unit
molecule = self.molecules[0]
charges = np.random.random([molecule.n_particles])
total_charge = molecule.total_charge
if type(total_charge) is unit.Quantity:
# Handle openforcefield >= 0.7.0
total_charge /= unit.elementary_charge
charges += (molecule.total_charge - charges.sum()) / molecule.n_particles
molecule.partial_charges = unit.Quantity(charges, unit.elementary_charge)
assert not np.all(molecule.partial_charges / unit.elementary_charge == 0)
assert (molecule.partial_charges is not None) and not np.all(molecule.partial_charges / unit.elementary_charge == 0)
# Add the molecule
generator.add_molecules(molecule)
# Create the System
Expand Down Expand Up @@ -707,7 +711,7 @@ def test_INSTALLED_FORCEFIELDS(self):
assert 'openff-1.1.0' in SMIRNOFFTemplateGenerator.INSTALLED_FORCEFIELDS
assert 'smirnoff99Frosst-1.1.0' in SMIRNOFFTemplateGenerator.INSTALLED_FORCEFIELDS
assert 'openff_unconstrained-1.1.0' not in SMIRNOFFTemplateGenerator.INSTALLED_FORCEFIELDS

def test_energies(self):
"""Test potential energies match between openforcefield and OpenMM ForceField"""
# DEBUG
Expand Down Expand Up @@ -753,7 +757,8 @@ def test_partial_charges_are_none(self):
from openforcefield.topology import Molecule
molecule = Molecule.from_smiles('C=O')
molecule.generate_conformers(n_conformers=1)
molecule._partial_charges = None
#molecule._partial_charges = None
assert (molecule.partial_charges is None) or np.all(molecule.partial_charges / unit.elementary_charge == 0)
# Test all supported SMIRNOFF force fields
for small_molecule_forcefield in SMIRNOFFTemplateGenerator.INSTALLED_FORCEFIELDS:
print(f'Testing energies for {small_molecule_forcefield}...')
Expand Down