Skip to content

Commit

Permalink
Check for enpty coordinates in RDKit converter (#4824)
Browse files Browse the repository at this point in the history
* prevents assignment of coordinates if all positions are zero
  • Loading branch information
RMeli authored Dec 9, 2024
1 parent 83f702f commit a27591c
Show file tree
Hide file tree
Showing 3 changed files with 22 additions and 6 deletions.
5 changes: 3 additions & 2 deletions package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -14,14 +14,15 @@ The rules for this file:


-------------------------------------------------------------------------------
??/??/?? IAlibay, ChiahsinChu
??/??/?? IAlibay, ChiahsinChu, RMeli

* 2.9.0

Fixes

Enhancements
* Added `precision` for XYZWriter (Issue #4775, PR #4771)
* Add check and warning for empty (all zero) coordinates in RDKit converter (PR #4824)
* Added `precision` for XYZWriter (Issue #4775, PR #4771)

Changes

Expand Down
9 changes: 6 additions & 3 deletions package/MDAnalysis/converters/RDKit.py
Original file line number Diff line number Diff line change
Expand Up @@ -363,9 +363,12 @@ def convert(self, obj, cache=True, NoImplicit=True, max_iter=200,

# add a conformer for the current Timestep
if hasattr(ag, "positions"):
if np.isnan(ag.positions).any():
warnings.warn("NaN detected in coordinates, the output "
"molecule will not have 3D coordinates assigned")
if np.isnan(ag.positions).any() or np.allclose(
ag.positions, 0.0, rtol=0.0, atol=1e-12
):
warnings.warn("NaN or empty coordinates detected in coordinates, "
"the output molecule will not have 3D coordinates "
"assigned")
else:
# assign coordinates
conf = Chem.Conformer(mol.GetNumAtoms())
Expand Down
14 changes: 13 additions & 1 deletion testsuite/MDAnalysisTests/converters/test_rdkit.py
Original file line number Diff line number Diff line change
Expand Up @@ -331,7 +331,7 @@ def test_nan_coords(self):
xyz = u.atoms.positions
xyz[0][2] = np.nan
u.atoms.positions = xyz
with pytest.warns(UserWarning, match="NaN detected"):
with pytest.warns(UserWarning, match="NaN .* detected"):
mol = u.atoms.convert_to("RDKIT")
with pytest.raises(ValueError, match="Bad Conformer Id"):
mol.GetConformer()
Expand Down Expand Up @@ -692,6 +692,18 @@ def test_reorder_atoms(self, smi):
expected = [a.GetSymbol() for a in mol.GetAtoms()]
assert values == expected

@pytest.mark.parametrize("smi", [
"O=S(C)(C)=NC",
])
def test_warn_empty_coords(self, smi):
mol = Chem.MolFromSmiles(smi)
mol = Chem.AddHs(mol)
# remove bond order and charges info
pdb = Chem.MolToPDBBlock(mol)
u = mda.Universe(StringIO(pdb), format="PDB")
with pytest.warns(match="NaN or empty coordinates detected"):
u.atoms.convert_to.rdkit()

def test_pdb_names(self):
u = mda.Universe(PDB_helix)
mol = u.atoms.convert_to.rdkit()
Expand Down

0 comments on commit a27591c

Please sign in to comment.