Skip to content

Commit

Permalink
First draft of new style Reader class
Browse files Browse the repository at this point in the history
Uses decorator to handle file handles

Addresses Issue #239
  • Loading branch information
richardjgowers committed Jun 14, 2015
1 parent 760eeb4 commit 0414ca1
Show file tree
Hide file tree
Showing 2 changed files with 67 additions and 33 deletions.
57 changes: 24 additions & 33 deletions package/MDAnalysis/coordinates/DLPoly.py
Original file line number Diff line number Diff line change
Expand Up @@ -131,7 +131,7 @@ def _read_first_frame(self):
ts._unitcell[2][:] = cellz


class HistoryReader(base.Reader):
class HistoryReader(base.MultiframeReader):
"""Reads DLPoly format HISTORY files
.. versionadded:: 0.11.0
Expand All @@ -154,49 +154,51 @@ def __init__(self, filename, convert_units=None, **kwargs):
self._dt = None
self._skip_timestep = None

# "private" file handle
self._file = open(self.filename, 'r')
self.title = self._file.readline().strip()
self._levcfg, self._imcon, self.numatoms = map(int, self._file.readline().split()[:3])
with open(self.filename, 'r') as f:
self.title = f.readline().strip()
self._levcfg, self._imcon, self.numatoms = map(int, f.readline().split()[:3])
self._startpoint = self._lastpos = f.tell()

self._has_vels = True if self._levcfg > 0 else False
self._has_forces = True if self._levcfg == 2 else False

self.ts = self._Timestep(self.numatoms,
velocities=self._has_vels,
forces=self._has_forces)
self._read_next_timestep()
self.next()

def _read_next_timestep(self, ts=None):
if ts is None:
ts = self.ts

line = self._file.readline() # timestep line
f = self._f

line = f.readline() # timestep line
if not line.startswith('timestep'):
raise IOError
if not self._imcon == 0:
ts._unitcell[0] = map(float, self._file.readline().split())
ts._unitcell[1] = map(float, self._file.readline().split())
ts._unitcell[2] = map(float, self._file.readline().split())
ts._unitcell[0] = map(float, f.readline().split())
ts._unitcell[1] = map(float, f.readline().split())
ts._unitcell[2] = map(float, f.readline().split())

# If ids are given, put them in here
# and later sort by them
ids = []

for i in range(self.numatoms):
line = self._file.readline().strip() # atom info line
line = f.readline().strip() # atom info line
try:
idx = int(line.split()[1])
except IndexError:
pass
else:
ids.append(idx)

# Read in this order for now, then later reorder in place
ts._pos[i] = map(float, self._file.readline().split())
ts._pos[i] = map(float, f.readline().split())
if self._has_vels:
ts._velocities[i] = map(float, self._file.readline().split())
ts._velocities[i] = map(float, f.readline().split())
if self._has_forces:
ts._forces[i] = map(float, self._file.readline().split())
ts._forces[i] = map(float, f.readline().split())

if ids:
ids = np.array(ids)
Expand All @@ -214,9 +216,9 @@ def _read_next_timestep(self, ts=None):

def _read_frame(self, frame):
"""frame is 0 based, error checking is done in base.getitem"""
self._file.seek(self._offsets[frame])
self._lastpos = self._offsets[frame]
self.ts.frame = frame # gets +1'd in read_next_frame
return self._read_next_timestep()
return self.next()

@property
def numframes(self):
Expand Down Expand Up @@ -260,25 +262,14 @@ def _read_numframes(self):

return numframes

def __iter__(self):
self._reopen()
while True:
try:
yield self._read_next_timestep()
except IOError:
self.rewind()
raise StopIteration

def rewind(self):
self._reopen()
self.next()

def _reopen(self):
self.close()
self._file = open(self.filename, 'r')
self._file.readline() # header is 2 lines
self._file.readline()
"""Change lastpos to just after the header,
this makes from_last_position decorator think the file has just
been opened
"""
self._lastpos = self._startpoint
self.ts.frame = 0

def close(self):
self._file.close()
43 changes: 43 additions & 0 deletions package/MDAnalysis/coordinates/base.py
Original file line number Diff line number Diff line change
Expand Up @@ -109,6 +109,7 @@
import warnings
import bisect
import numpy as np
import functools

from MDAnalysis.core import units, flags
from MDAnalysis.core.util import asiterable
Expand Down Expand Up @@ -1195,3 +1196,45 @@ def close(self):

def __del__(self):
self.close()


def from_last_position(func):
@functools.wraps(func)
def wrapper(self, *args, **kwargs):
with open(self.filename, 'r') as self._f:
self._f.seek(self._lastpos)
ret = func(self, *args, **kwargs)
self._lastpos = self._f.tell()
return ret
return wrapper

class MultiframeReader(Reader):
"""Designed to not require an open file handle at all times
@from_last_position decorator opens the file handle at the last known
position and saves the last known position when done
"""
@from_last_position
def next(self):
return self._read_next_timestep()

def __iter__(self):
self._reopen()
while True:
try:
yield self.next()
except IOError:
self.rewind()
raise StopIteration

def close(self):
"""There are never any file handles left open to close"""
pass

def __getstate__(self):
result = self.__dict__.copy()
del result['_f'] # closed file handle
return result

def __setstate__(self, d):
self.__dict__ = d

0 comments on commit 0414ca1

Please sign in to comment.