Skip to content

Commit

Permalink
adding element attribute to txyz parser (#3826)
Browse files Browse the repository at this point in the history
* adding element attribute to txyz parser

* adding elements when names is valid only

* update docstring and change log

* Update TXYZParser.py

* Update TXYZParser.py

* Update TXYZParser.py

* Update TXYZParser.py

* Update TXYZParser.py

* Update TXYZParser.py

* Update TXYZParser.py

* Update CHANGELOG

* Apply suggestions from code review

Co-authored-by: Jonathan Barnoud <[email protected]>

* Apply suggestions from code review

Co-authored-by: Jonathan Barnoud <[email protected]>

Co-authored-by: Jonathan Barnoud <[email protected]>
  • Loading branch information
aya9aladdin and jbarnoud authored Sep 21, 2022
1 parent 80ff647 commit fa8b03b
Show file tree
Hide file tree
Showing 3 changed files with 42 additions and 5 deletions.
2 changes: 1 addition & 1 deletion package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -13,7 +13,6 @@ The rules for this file:
* release numbers follow "Semantic Versioning" http://semver.org

------------------------------------------------------------------------------

??/??/?? IAlibay, Luthaf, hmacdope, rafaelpap, jbarnoud, BFedder, aya9aladdin

* 2.4.0
Expand Down Expand Up @@ -53,6 +52,7 @@ Changes
* The "AMBER" entry in setup.py's extra_requires has now been removed. Please
use `pip install ./package[extra_formats]` instead (PR #3810)
* `Universe.empty` emmits less warnings (PR #3814)
* adding element attribute to TXYZParser if all atom names are valid element symbols (PR #3826)

Deprecations

Expand Down
15 changes: 14 additions & 1 deletion package/MDAnalysis/topology/TXYZParser.py
Original file line number Diff line number Diff line change
Expand Up @@ -45,8 +45,10 @@

import itertools
import numpy as np
import warnings

from . import guessers
from .tables import SYMB2Z
from ..lib.util import openany
from .base import TopologyReaderBase
from ..core.topology import Topology
Expand All @@ -59,6 +61,7 @@
Resids,
Resnums,
Segids,
Elements,
)


Expand All @@ -69,8 +72,11 @@ class TXYZParser(TopologyReaderBase):
- Atomnames
- Atomtypes
- Elements (if all atom names are element symbols)
.. versionadded:: 0.17.0
.. versionchanged:: 2.4.0
Adding the `Element` attribute if all names are valid element symbols.
"""
format = ['TXYZ', 'ARC']

Expand Down Expand Up @@ -126,7 +132,14 @@ def parse(self, **kwargs):
Resnums(np.array([1])),
Segids(np.array(['SYSTEM'], dtype=object)),
]

if all(n.capitalize() in SYMB2Z for n in names):
attrs.append(Elements(np.array(names, dtype=object)))

else:
warnings.warn("Element information is missing, elements attribute "
"will not be populated. If needed these can be "
"guessed using MDAnalysis.topology.guessers.")

top = Topology(natoms, 1, 1,
attrs=attrs)

Expand Down
30 changes: 27 additions & 3 deletions testsuite/MDAnalysisTests/topology/test_txyz.py
Original file line number Diff line number Diff line change
Expand Up @@ -22,15 +22,16 @@
#
import MDAnalysis as mda
import pytest
import numpy as np

from MDAnalysisTests.topology.base import ParserBase
from MDAnalysisTests.datafiles import TXYZ, ARC

from MDAnalysisTests.datafiles import TXYZ, ARC, ARC_PBC
from numpy.testing import assert_equal

class TestTXYZParser(ParserBase):
parser = mda.topology.TXYZParser.TXYZParser
guessed_attrs = ['masses']
expected_attrs = ['ids', 'names', 'bonds', 'types']
expected_attrs = ['ids', 'names', 'bonds', 'types', 'elements']
expected_n_residues = 1
expected_n_atoms = 9
expected_n_segments = 1
Expand All @@ -51,3 +52,26 @@ def test_atom_type_type(self, top):
types = top.types.values
type_is_str = [isinstance(atom_type, str) for atom_type in types]
assert all(type_is_str)

def test_TXYZ_elements():
"""The test checks whether elements attribute are assigned
properly given a TXYZ file with valid elements record.
"""
u = mda.Universe(TXYZ, format='TXYZ')
element_list = np.array(['C', 'H', 'H', 'O', 'H', 'C', 'H', 'H', 'H'], dtype=object)
assert_equal(u.atoms.elements, element_list)


def test_missing_elements_noattribute():
"""Check that:
1) a warning is raised if elements are missing
2) the elements attribute is not set
"""
wmsg = ("Element information is missing, elements attribute will not be "
"populated")
with pytest.warns(UserWarning, match=wmsg):
u = mda.Universe(ARC_PBC)
with pytest.raises(AttributeError):
_ = u.atoms.elements

0 comments on commit fa8b03b

Please sign in to comment.