Skip to content

Commit

Permalink
Merge branch 'develop' into distopia_0.3.0_compat
Browse files Browse the repository at this point in the history
  • Loading branch information
orbeckst authored Jan 6, 2025
2 parents 730cee9 + 536a390 commit a865a46
Show file tree
Hide file tree
Showing 5 changed files with 76 additions and 1 deletion.
1 change: 1 addition & 0 deletions package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -29,6 +29,7 @@ Fixes
Enhancements
* Improve distopia backend support in line with new functionality available
in distopia >= 0.3.1 (PR #4734)
* Addition of 'water' token for water selection (Issue #4839)
* Enables parallelization for analysis.density.DensityAnalysis (Issue #4677, PR #4729)
* Enables parallelization for analysis.contacts.Contacts (Issue #4660)
* Enable parallelization for analysis.nucleicacids.NucPairDist (Issue #4670)
Expand Down
31 changes: 31 additions & 0 deletions package/MDAnalysis/core/selection.py
Original file line number Diff line number Diff line change
Expand Up @@ -1094,6 +1094,37 @@ def _apply(self, group):
return group[mask]


class WaterSelection(Selection):
"""All atoms in water residues with recognized water residue names.
Recognized residue names:
* recognized 3 Letter resnames: 'H2O', 'HOH', 'OH2', 'HHO', 'OHH'
'TIP', 'T3P', 'T4P', 'T5P', 'SOL', 'WAT'
* recognized 4 Letter resnames: 'TIP2', 'TIP3', 'TIP4'
.. versionadded:: 2.9.0
"""
token = 'water'

# Recognized water resnames
water_res = {
'H2O', 'HOH', 'OH2', 'HHO', 'OHH',
'T3P', 'T4P', 'T5P', 'SOL', 'WAT',
'TIP', 'TIP2', 'TIP3', 'TIP4'
}

def _apply(self, group):
resnames = group.universe._topology.resnames
nmidx = resnames.nmidx[group.resindices]

matches = [ix for (nm, ix) in resnames.namedict.items()
if nm in self.water_res]
mask = np.isin(nmidx, matches)

return group[mask]

class BackboneSelection(ProteinSelection):
"""A BackboneSelection contains all atoms with name 'N', 'CA', 'C', 'O'.
Expand Down
Original file line number Diff line number Diff line change
@@ -1,3 +1,3 @@
.. automodule:: MDAnalysis.analysis.polymer
:members:

:inherited-members:
9 changes: 9 additions & 0 deletions package/doc/sphinx/source/documentation_pages/selections.rst
Original file line number Diff line number Diff line change
Expand Up @@ -114,6 +114,15 @@ protein, backbone, nucleic, nucleicbackbone
is identfied by a hard-coded set of residue names so it may not
work for esoteric residues.

water
selects all atoms that belong to a set of water residues; water
is defined with a set of common water abbreviations present in
topology files and may not work with certain water residue names.
Currently the following water resnames are supported:
3 letter resnames: ``H2O``, ``HOH``, ``OH2``, ``HHO``, ``OHH``, ``TIP``,
``T3P``, ``T4P``, ``T5P``, ``SOL``, ``WAT``.
4 letter resnames: ``TIP2``, ``TIP3``, ``TIP4``.

segid *seg-name*
select by segid (as given in the topology), e.g. ``segid 4AKE`` or
``segid DMPC``
Expand Down
34 changes: 34 additions & 0 deletions testsuite/MDAnalysisTests/core/test_atomselections.py
Original file line number Diff line number Diff line change
Expand Up @@ -48,6 +48,8 @@
PDB_helix,
PDB_elements,
PDB_charges,
waterPSF,
PDB_full,
)
from MDAnalysisTests import make_Universe

Expand Down Expand Up @@ -650,6 +652,38 @@ def test_nucleicsugar(self, universe):
assert_equal(rna.n_atoms, rna.n_residues * 5)


class TestSelectionsWater(object):
@pytest.fixture(scope='class')
def universe(self):
return MDAnalysis.Universe(GRO)

@pytest.fixture(scope='class')
def universe2(self):
return MDAnalysis.Universe(waterPSF)

@pytest.fixture(scope='class')
def universe3(self):
return MDAnalysis.Universe(PDB_full)

def test_water_gro(self, universe):
# Test SOL water with 4 atoms
water_gro = universe.select_atoms("water")
assert_equal(water_gro.n_atoms, 44336)
assert_equal(water_gro.n_residues, 11084)

def test_water_tip3(self, universe2):
# Test TIP3 water with 3 atoms
water_tip3 = universe2.select_atoms('water')
assert_equal(water_tip3.n_atoms, 15)
assert_equal(water_tip3.n_residues, 5)

def test_water_pdb(self, universe3):
# Test HOH water with 1 atom
water_pdb = universe3.select_atoms("water")
assert_equal(water_pdb.n_residues, 188)
assert_equal(water_pdb.n_atoms, 188)


class BaseDistanceSelection(object):
"""Both KDTree and distmat selections on orthogonal system
Expand Down

0 comments on commit a865a46

Please sign in to comment.