Skip to content

Commit

Permalink
fix error message for XDR readers
Browse files Browse the repository at this point in the history
- error message for not being able to write lockfile did not insert the filename into
  message: fixed by using proper f string
- changed all .format to modern f string in XDR.py
- add explicit test for offsets_filename()
- add additional pytest.warns() to tests to catch warnings (reduces warnings output and
  checks that all expected warnings are raised with proper messages)
  • Loading branch information
orbeckst committed Aug 10, 2023
1 parent d33cc3a commit dc10474
Show file tree
Hide file tree
Showing 2 changed files with 35 additions and 17 deletions.
17 changes: 8 additions & 9 deletions package/MDAnalysis/coordinates/XDR.py
Original file line number Diff line number Diff line change
Expand Up @@ -63,8 +63,7 @@ def offsets_filename(filename, ending='npz'):
"""
head, tail = split(filename)
return join(head, '.{tail}_offsets.{ending}'.format(tail=tail,
ending=ending))
return join(head, f'.{tail}_offsets.{ending}')


def read_numpy_offsets(filename):
Expand All @@ -88,7 +87,7 @@ def read_numpy_offsets(filename):

# `ValueError` is encountered when the offset file is corrupted.
except (ValueError, IOError):
warnings.warn("Failed to load offsets file {}\n".format(filename))
warnings.warn(f"Failed to load offsets file {filename}\n")
return False

class XDRBaseReader(base.ReaderBase):
Expand Down Expand Up @@ -201,7 +200,7 @@ def _load_offsets(self):
except OSError as e:
if isinstance(e, PermissionError) or e.errno == errno.EROFS:
warnings.warn(f"Cannot write lock/offset file in same location as "
"{self.filename}. Using slow offset calculation.")
f"{self.filename}. Using slow offset calculation.")
self._read_offsets(store=True)
return
else:
Expand All @@ -220,10 +219,10 @@ def _load_offsets(self):
# refer to Issue #1893
data = read_numpy_offsets(fname)
if not data:
warnings.warn("Reading offsets from {} failed, "
"reading offsets from trajectory instead.\n"
"Consider setting 'refresh_offsets=True' "
"when loading your Universe.".format(fname))
warnings.warn(f"Reading offsets from {fname} failed, "
f"reading offsets from trajectory instead.\n"
f"Consider setting 'refresh_offsets=True' "
f"when loading your Universe.")
self._read_offsets(store=True)
return

Expand Down Expand Up @@ -256,7 +255,7 @@ def _read_offsets(self, store=False):
offsets=offsets, size=size, ctime=ctime,
n_atoms=self._xdr.n_atoms)
except Exception as e:
warnings.warn("Couldn't save offsets because: {}".format(e))
warnings.warn(f"Couldn't save offsets because: {e}")

@property
def n_frames(self):
Expand Down
35 changes: 27 additions & 8 deletions testsuite/MDAnalysisTests/coordinates/test_xdr.py
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,7 @@
import pytest
from unittest.mock import patch

import re
import os
import shutil
import subprocess
Expand All @@ -46,6 +47,14 @@
from MDAnalysis.coordinates import XDR
from MDAnalysisTests.util import get_userid

@pytest.mark.parametrize("filename,kwargs,reference", [
("foo.xtc", {}, ".foo.xtc_offsets.npz"),
("foo.xtc", {"ending": "npz"}, ".foo.xtc_offsets.npz"),
("bar.0001.trr", {"ending": "npzzzz"}, ".bar.0001.trr_offsets.npzzzz"),
])
def test_offsets_filename(filename, kwargs, reference):
fn = XDR.offsets_filename(filename, **kwargs)
assert fn == reference

class _XDRReader_Sub(object):
@pytest.fixture()
Expand Down Expand Up @@ -212,7 +221,7 @@ def go_beyond_EOF():

with pytest.raises(StopIteration):
go_beyond_EOF()

def test_read_next_timestep_ts_no_positions(self, universe):
# primarily tests branching on ts.has_positions in _read_next_timestep
ts = universe.trajectory[0]
Expand Down Expand Up @@ -755,25 +764,34 @@ def test_nonexistent_offsets_file(self, traj):
outfile_offsets = XDR.offsets_filename(traj)
with patch.object(np, "load") as np_load_mock:
np_load_mock.side_effect = IOError
saved_offsets = XDR.read_numpy_offsets(outfile_offsets)
assert_equal(saved_offsets, False)
with pytest.warns(UserWarning, match=re.escape(
f"Failed to load offsets file {outfile_offsets}")):
saved_offsets = XDR.read_numpy_offsets(outfile_offsets)
assert saved_offsets == False

def test_nonexistent_offsets_file(self, traj):
def test_corrupted_offsets_file(self, traj):
# assert that a corrupted file returns False during read-in
# Issue #3230
outfile_offsets = XDR.offsets_filename(traj)
with patch.object(np, "load") as np_load_mock:
np_load_mock.side_effect = ValueError
saved_offsets = XDR.read_numpy_offsets(outfile_offsets)
assert_equal(saved_offsets, False)
with pytest.warns(UserWarning, match=re.escape(
f"Failed to load offsets file {outfile_offsets}")):
saved_offsets = XDR.read_numpy_offsets(outfile_offsets)
assert saved_offsets == False

def test_reload_offsets_if_offsets_readin_io_fails(self, trajectory):
# force the np.load call that is called in read_numpy_offsets
# during _load_offsets to give an IOError
# ensure that offsets are then read-in from the trajectory
with patch.object(np, "load") as np_load_mock:
np_load_mock.side_effect = IOError
trajectory._load_offsets()
with pytest.warns(UserWarning) as record:
trajectory._load_offsets()
assert len(record) == 2
assert str(record[0].message).startswith("Failed to load offsets file")
assert "reading offsets from trajectory instead" in str(record[1].message)

assert_almost_equal(
trajectory._xdr.offsets,
self.ref_offsets,
Expand Down Expand Up @@ -866,7 +884,8 @@ def test_unsupported_format(self, traj):
np.savez(fname, **saved_offsets)

# ok as long as this doesn't throw
reader = self._reader(traj)
with pytest.warns(UserWarning, match="Reload offsets from trajectory"):
reader = self._reader(traj)
reader[idx_frame]

@pytest.mark.skipif(get_userid() == 0, reason="cannot readonly as root")
Expand Down

0 comments on commit dc10474

Please sign in to comment.