Skip to content

Commit

Permalink
warnings for ALL instant selectors
Browse files Browse the repository at this point in the history
Thanks to @richardjgowers for insights #1403 (comment)
  • Loading branch information
orbeckst committed Jun 17, 2017
1 parent 75eade4 commit eca6a61
Show file tree
Hide file tree
Showing 3 changed files with 74 additions and 19 deletions.
33 changes: 15 additions & 18 deletions package/MDAnalysis/core/groups.py
Original file line number Diff line number Diff line change
Expand Up @@ -1504,9 +1504,7 @@ class AtomGroup(GroupBase):
.. deprecated:: 0.16.2
*Instant selectors* of AtomGroup will be removed in the 1.0 release.
See issue `#1377
<https://github.com/MDAnalysis/mdanalysis/issues/1377>`_ for
more details.
See :ref:`Instant selectors`_ for details and alternatives.
"""
def __getitem__(self, item):
Expand All @@ -1515,10 +1513,6 @@ def __getitem__(self, item):
#
# u.atoms['HT1'] access, otherwise default
if isinstance(item, string_types):
warnings.warn("Instant selectors AtomGroup['<name>'] "
"will be removed in 1.0. "
"Use `select_atoms('name <name>')`.",
DeprecationWarning)
try:
return self._get_named_atom(item)
except (AttributeError, selection.SelectionError):
Expand All @@ -1536,9 +1530,6 @@ def __getattr__(self, attr):
# raise NDE(_ATTR_ERRORS[attr])
raise NoDataError("AtomGroup has no fragments; this requires Bonds")
elif hasattr(self.universe._topology, 'names'):
warnings.warn("Instant selectors AtomGroup.<name> "
"will be removed in 1.0",
DeprecationWarning)
# Ugly hack to make multiple __getattr__s work
try:
return self._get_named_atom(attr)
Expand Down Expand Up @@ -2208,6 +2199,11 @@ class ResidueGroup(GroupBase):
ResidueGroups can be compared and combined using group operators. See the
list of these operators on :class:`GroupBase`.
.. deprecated:: 0.16.2
*Instant selectors* of Segments will be removed in the 1.0 release.
See :ref:`Instant selectors`_ for details and alternatives.
"""

@property
Expand Down Expand Up @@ -2322,6 +2318,11 @@ class SegmentGroup(GroupBase):
SegmentGroups can be compared and combined using group operators. See the
list of these operators on :class:`GroupBase`.
.. deprecated:: 0.16.2
*Instant selectors* of Segments will be removed in the 1.0 release.
See :ref:`Instant selectors`_ for details and alternatives.
"""

@property
Expand Down Expand Up @@ -2669,9 +2670,7 @@ class Segment(ComponentBase):
.. deprecated:: 0.16.2
*Instant selectors* of Segments will be removed in the 1.0 release.
See issue `#1377
<https://github.com/MDAnalysis/mdanalysis/issues/1377>`_ for
more details.
See :ref:`Instant selectors`_ for details and alternatives.
"""
def __repr__(self):
Expand All @@ -2696,16 +2695,14 @@ def __getattr__(self, attr):
if attr.startswith('r') and attr[1:].isdigit():
resnum = int(attr[1:])
rg = self.residues[resnum - 1] # convert to 0 based
warnings.warn("Instant selectors Segment.r<resid> will be removed in 1.0",
warnings.warn("Instant selectors Segment.r<N> will be removed in 1.0. "
"Use Segment.residues[N-1] instead.",
DeprecationWarning)
return rg
# Resname accesss
if hasattr(self.residues, 'resnames'):
try:
rg = self.residues._get_named_residue(attr)
warnings.warn("Instant selectors Segment.<resname> will be removed in 1.0",
DeprecationWarning)
return rg
return self.residues._get_named_residue(attr)
except selection.SelectionError:
pass
raise AttributeError("{cls} has no attribute {attr}"
Expand Down
40 changes: 40 additions & 0 deletions package/MDAnalysis/core/topologyattrs.py
Original file line number Diff line number Diff line change
Expand Up @@ -42,6 +42,8 @@
import numbers
import numpy as np

from numpy.lib.utils import deprecate

from . import flags
from ..lib.util import cached, convert_aa_code, iterable
from ..lib import transformations, mdamath
Expand Down Expand Up @@ -409,6 +411,9 @@ def getattr__(atomgroup, name):
raise AttributeError("'{0}' object has no attribute '{1}'".format(
atomgroup.__class__.__name__, name))

@deprecate(message="Instant selector AtomGroup['<name>'] or AtomGroup.<name> "
"is deprecated and will be removed in 1.0. "
"Use AtomGroup.select_atoms('name <name>') instead.")
def _get_named_atom(group, name):
"""Get all atoms with name *name* in the current AtomGroup.
Expand All @@ -417,6 +422,14 @@ def _get_named_atom(group, name):
no atoms are found, a :exc:`SelectionError` is raised.
.. versionadded:: 0.9.2
.. deprecated:: 0.16.2
*Instant selectors* will be removed in the 1.0 release.
Use ``AtomGroup.select_atoms('name <name>')`` instead.
See issue `#1377
<https://github.com/MDAnalysis/mdanalysis/issues/1377>`_ for
more details.
"""
# There can be more than one atom with the same name
atomlist = group.atoms.unique[group.atoms.unique.names == name]
Expand Down Expand Up @@ -1057,6 +1070,13 @@ def getattr__(residuegroup, resname):
# This transplant is hardcoded for now to allow for multiple getattr things
#transplants[Segment].append(('__getattr__', getattr__))


@deprecate(message="Instant selector ResidueGroup.<name> "
"or Segment.<name> "
"is deprecated and will be removed in 1.0. "
"Use ResidueGroup[ResidueGroup.resnames == '<name>'] "
"or Segment.residues[Segment.residues == '<name>'] "
"instead.")
def _get_named_residue(group, resname):
"""Get all residues with name *resname* in the current ResidueGroup
or Segment.
Expand All @@ -1068,6 +1088,15 @@ def _get_named_residue(group, resname):
.. versionadded:: 0.9.2
.. deprecated:: 0.16.2
*Instant selectors* will be removed in the 1.0 release.
Use ``ResidueGroup[ResidueGroup.resnames == '<name>']``
or ``Segment.residues[Segment.residues == '<name>']``
instead.
See issue `#1377
<https://github.com/MDAnalysis/mdanalysis/issues/1377>`_ for
more details.
"""
# There can be more than one residue with the same name
residues = group.residues.unique[
Expand Down Expand Up @@ -1244,6 +1273,10 @@ def getattr__(segmentgroup, segid):
transplants[SegmentGroup].append(
('__getattr__', getattr__))

@deprecate(message="Instant selector SegmentGroup.<name> "
"is deprecated and will be removed in 1.0. "
"Use SegmentGroup[SegmentGroup.segids == '<name>'] "
"instead.")
def _get_named_segment(group, segid):
"""Get all segments with name *segid* in the current SegmentGroup.
Expand All @@ -1254,6 +1287,13 @@ def _get_named_segment(group, segid):
.. versionadded:: 0.9.2
.. deprecated:: 0.16.2
*Instant selectors* will be removed in the 1.0 release.
Use ``SegmentGroup[SegmentGroup.segids == '<name>']`` instead.
See issue `#1377
<https://github.com/MDAnalysis/mdanalysis/issues/1377>`_ for
more details.
"""
# Undo adding 's' if segid started with digit
if segid.startswith('s') and len(segid) >= 2 and segid[1].isdigit():
Expand Down
20 changes: 19 additions & 1 deletion package/doc/sphinx/source/documentation_pages/selections.rst
Original file line number Diff line number Diff line change
Expand Up @@ -286,7 +286,6 @@ Instant selectors
=================

.. deprecated:: 0.16.2

*Instant selectors* will be removed in the 1.0 release in order to
streamline the MDAnalysis user interface. They do not seem to be
widely used anymore, can produce cryptic error messages, and are
Expand All @@ -309,6 +308,11 @@ other levels of the structural hierarchy, namely for
Segment selector
----------------

.. deprecated:: 0.16.2
Use ``SegmentGroup[SegmentGroup.segids == '<name>']`` instead. Note that this
*always* returns a :class:`SegmentGroup` and *never* a :class:`Segment`
(unlike the instant selector).

- ``universe.<segid>`` or ``universe.s<segid>`` (if *<segid>* starts with a
number)
- returns a :class:`~MDAnalysis.core.groups.Segment`
Expand All @@ -320,6 +324,9 @@ Segment selector
Resid selector
--------------

.. deprecated:: 0.16.2
Use ``Segment.residues[N-1]`` instead.

- ``seg.r<N>`` selects residue with number ``<N>``
- returns a :class:`~MDAnalysis.core.groups.Residue`
- works for :class:`~MDAnalysis.core.groups.Segment` and :class:`~MDAnalysis.core.groups.SegmentGroup`
Expand All @@ -330,6 +337,12 @@ Resid selector
Residue name selector
---------------------

.. deprecated:: 0.16.2
Use ``ResidueGroup[ResidueGroup.resnames == '<name>']`` or
``Segment.residues[Segment.residues == '<name>']`` instead. Note that this
*always* returns a :class:`ResidueGroup` and *never* a :class:`Residue`
(unlike the instant selector).

- ``seg.<resname>`` selects residues with residue name ``<resname>``
- returns a :class:`~MDAnalysis.core.groups.ResidueGroup`
- works for :class:`~MDAnalysis.core.groups.Segment` and :class:`~MDAnalysis.core.groups.SegmentGroup`
Expand All @@ -346,6 +359,11 @@ Residue name selector
Atom name selector
------------------

.. deprecated:: 0.16.2
Use ``AtomGroup.select_atoms('name <name>')`` instead. Note that this
*always* returns an :class:`AtomGroup` and *never* an :class:`Atom` (unlike
the instant selector).

- ``g.<atomname>`` selects a single atom or a group of atoms with name
``<atomname>``
- returns
Expand Down

0 comments on commit eca6a61

Please sign in to comment.