Skip to content

Commit

Permalink
testing flavor of DCDReader/Writer
Browse files Browse the repository at this point in the history
Introduced an experimental class-attribute flavor for the DCDReader/DCDWriter that
helps to distinguish the CHARMM DCD code from the LAMMPS code (and later will be
useful for the NAMD DCD Reader/Writer in #359 )

Added tests to check that the LAMMPS DCDWriter and Reader are being used because
the last coverage report said that the LAMMPS.DCDWriter was not being covered.
  • Loading branch information
orbeckst committed Oct 6, 2015
1 parent 15d577b commit f5ad35f
Show file tree
Hide file tree
Showing 3 changed files with 72 additions and 6 deletions.
2 changes: 2 additions & 0 deletions package/MDAnalysis/coordinates/DCD.py
Original file line number Diff line number Diff line change
Expand Up @@ -192,6 +192,7 @@ class DCDWriter(base.Writer):
.. _Issue 187: https://github.com/MDAnalysis/mdanalysis/issues/187
"""
format = 'DCD'
flavor = 'CHARMM'
units = {'time': 'AKMA', 'length': 'Angstrom'}

def __init__(self, filename, n_atoms, start=0, step=1,
Expand Down Expand Up @@ -399,6 +400,7 @@ class DCDReader(base.Reader):
Removed skip keyword and functionality
"""
format = 'DCD'
flavor = 'CHARMM'
units = {'time': 'AKMA', 'length': 'Angstrom'}
_Timestep = Timestep

Expand Down
8 changes: 5 additions & 3 deletions package/MDAnalysis/coordinates/LAMMPS.py
Original file line number Diff line number Diff line change
Expand Up @@ -79,10 +79,11 @@ class DCDWriter(DCD.DCDWriter):
"Angstrom". See :mod:`MDAnalysis.units` for other recognized
values.
"""
format = "DCD"
format = 'DCD'
flavor = 'LAMMPS'

def __init__(self, *args, **kwargs):
self.units = {'time': 'ps', 'length': 'Angstrom'} # must be instance level
self.units = {'time': 'fs', 'length': 'Angstrom'} # must be instance level
self.units['time'] = kwargs.pop('timeunit', self.units['time'])
self.units['length'] = kwargs.pop('lengthunit', self.units['length'])
for unit_type, unit in self.units.items():
Expand All @@ -105,7 +106,8 @@ class DCDReader(DCD.DCDReader):
.. _units style: http://lammps.sandia.gov/doc/units.html
..
"""
format = "DCD"
format = 'DCD'
flavor = 'LAMMPS'

def __init__(self, dcdfilename, **kwargs):
self.units = {'time': 'fs', 'length': 'Angstrom'} # must be instance level
Expand Down
68 changes: 65 additions & 3 deletions testsuite/MDAnalysisTests/coordinates/test_lammps.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,8 @@

import MDAnalysis as mda

from numpy.testing import (assert_equal, assert_almost_equal, assert_raises)
from numpy.testing import (assert_equal, assert_almost_equal, assert_raises,
assert_)
import tempdir
from unittest import TestCase

Expand All @@ -12,7 +13,7 @@
RefLAMMPSDataDCD)


def test_datareader_VE():
def test_datareader_ValueError():
from MDAnalysis.coordinates.LAMMPS import DATAReader
assert_raises(ValueError, DATAReader, 'filename')

Expand Down Expand Up @@ -66,13 +67,18 @@ class TestLammpsDataMini_Coords(_TestLammpsData_Coords, RefLAMMPSDataMini):
# need more tests of the LAMMPS DCDReader

class TestLAMMPSDCDReader(TestCase, RefLAMMPSDataDCD):
flavor = 'LAMMPS'

def setUp(self):
self.u = mda.Universe(self.topology, self.trajectory,
format=self.format)

def tearDown(self):
del self.u

def test_Reader_is_LAMMPS(self):
assert_(self.u.trajectory.flavor, self.flavor)

def get_frame_from_end(self, offset):
iframe = self.u.trajectory.n_frames - 1 - offset
iframe = iframe if iframe > 0 else 0
Expand Down Expand Up @@ -123,6 +129,8 @@ def wrong_load(unit="GARBAGE"):


class TestLAMMPSDCDWriter(TestCase, RefLAMMPSDataDCD):
flavor = 'LAMMPS'

def setUp(self):
self.u = mda.Universe(self.topology, self.trajectory,
format=self.format)
Expand All @@ -139,6 +147,11 @@ def tearDown(self):
del self.u
del self.tmpdir

def test_Writer_is_LAMMPS(self):
with mda.Writer(self.outfile, n_atoms=self.u.atoms.n_atoms,
format=self.format) as W:
assert_(W.flavor, self.flavor)

def test_Writer(self, n_frames=3):
W = mda.Writer(self.outfile, n_atoms=self.u.atoms.n_atoms,
format=self.format)
Expand All @@ -152,7 +165,9 @@ def test_Writer(self, n_frames=3):
self.u.trajectory[n_frames - 1].positions,
6, err_msg="coordinate mismatch between corresponding frames")


def test_OtherWriter_is_LAMMPS(self):
with self.u.trajectory.OtherWriter(self.outfile) as W:
assert_(W.flavor, self.flavor)

def test_OtherWriter(self):
times = []
Expand All @@ -174,3 +189,50 @@ def test_OtherWriter(self):
assert_almost_equal(reversed.trajectory[-1].positions,
self.u.trajectory[0].positions,
6, err_msg="coordinate mismatch between corresponding frames")

class TestLAMMPSDCDWriterClass(TestCase):
flavor = 'LAMMPS'

def setUp(self):
# dummy output file
ext = ".dcd"
self.tmpdir = tempdir.TempDir()
self.outfile = os.path.join(self.tmpdir.name, 'lammps-writer-test' + ext)

def tearDown(self):
try:
os.unlink(self.outfile)
except:
pass
del self.tmpdir

def test_Writer_is_LAMMPS(self):
with mda.coordinates.LAMMPS.DCDWriter(self.outfile, n_atoms=10) as W:
assert_(W.flavor, self.flavor)

def test_open(self):
def open_dcd():
try:
with mda.coordinates.LAMMPS.DCDWriter(self.outfile, n_atoms=10):
pass
except Exception:
return False
else:
return True
assert_(open_dcd(), True)

def test_wrong_time_unit(self):
def wrong_load(unit="nm"):
with mda.coordinates.LAMMPS.DCDWriter(self.outfile, n_atoms=10,
timeunit=unit):
pass
assert_raises(TypeError, wrong_load)

def test_wrong_unit(self):
def wrong_load(unit="GARBAGE"):
with mda.coordinates.LAMMPS.DCDWriter(self.outfile, n_atoms=10,
timeunit=unit):
pass
assert_raises(ValueError, wrong_load)


2 comments on commit f5ad35f

@richardjgowers
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

https://github.com/MDAnalysis/mdanalysis/blob/develop/testsuite/MDAnalysisTests/test_util.py#L547

I should probably have added Writers to the format checking here

@orbeckst
Copy link
Member Author

@orbeckst orbeckst commented on f5ad35f Oct 6, 2015 via email

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Please sign in to comment.