Skip to content

Commit

Permalink
Updated Reader API
Browse files Browse the repository at this point in the history
fixed, skip_timestep, periodic, skip and delta no longer required for
Readers (Issue #350)

All writers now refer to time between steps as "dt" (was previously
delta in XTC, TRR and NCDF Writers) (Issue #206)

Timesteps now have a default dt of 1.0 ps (and issue warning when using
this)

Timestep attribute time now tries to return data['time'] first, then dt
* frame, using the default value of 1.0 ps if necessary. (Issue #308)
  • Loading branch information
richardjgowers committed Jul 19, 2015
1 parent 606454d commit 0e36461
Show file tree
Hide file tree
Showing 16 changed files with 176 additions and 201 deletions.
8 changes: 8 additions & 0 deletions package/CHANGELOG
Original file line number Diff line number Diff line change
Expand Up @@ -74,8 +74,16 @@ The rules for this file:
should not affect any code because
MDAnalysis.selections.get_writer() already existed.
* Moved norm, normal, angle, stp and dihedral from lib.util to lib.mdamath
* All Readers have a default "dt" of 1.0 ps. This applies also to single
frame readers, these would previously raise an error on accessing dt.
(Issue #350)
* NCDF Reader no longer has a default skip_timestep of 1 (Issue #350)
* All Writers now refer to time between written Timesteps as "dt",
was previously "delta" in some Writers. (Issue #206)

Fixes
* Amber TRJ and NCDF Reader & Writers now use 'dt' instead of 'delta'
to refer to time elapsed between timesteps. (Issue #350)
* Fixed TPRParser considering LJ 1..4 exclusions as bonds (Issue #351)
* ChainReaders no longer leak (Issue #312)
* analysis.hbonds.HydrogenBondAnalysis performs a sanity check for
Expand Down
12 changes: 11 additions & 1 deletion package/MDAnalysis/coordinates/DCD.py
Original file line number Diff line number Diff line change
Expand Up @@ -419,7 +419,13 @@ def __init__(self, dcdfilename, **kwargs):
self.fixed = 0
self.periodic = False

# This reads skip_timestep and delta from header
self._read_dcd_header()

# Convert delta to ps
delta = mdaunits.convert(self.delta, 'AKMA', 'ps')

self._ts_kwargs.update({'dt':self.skip_timestep * delta})
self.ts = self._Timestep(self.numatoms, **self._ts_kwargs)
# Read in the first timestep
self._read_next_timestep()
Expand Down Expand Up @@ -472,7 +478,7 @@ def _read_next_timestep(self, ts=None):
"""
if ts is None:
ts = self.ts
ts._frame = self._read_next_frame(ts._x, ts._y, ts._z, ts._unitcell, self.skip)
ts._frame = self._read_next_frame(ts._x, ts._y, ts._z, ts._unitcell, 1)

ts.frame += 1
return ts
Expand Down Expand Up @@ -583,6 +589,10 @@ def Writer(self, filename, **kwargs):
# dt keyword is simply passed through if provided
return DCDWriter(filename, numatoms, **kwargs)

@property
def dt(self):
return self.skip_timestep * self.convert_time_from_native(self.delta)


DCDReader._read_dcd_header = new.instancemethod(_dcdmodule.__read_dcd_header, None, DCDReader)
DCDReader._read_next_frame = new.instancemethod(_dcdmodule.__read_next_frame, None, DCDReader)
Expand Down
1 change: 0 additions & 1 deletion package/MDAnalysis/coordinates/DLPoly.py
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,6 @@

import numpy as np

import MDAnalysis
from . import base
from . import core

Expand Down
5 changes: 1 addition & 4 deletions package/MDAnalysis/coordinates/GMS.py
Original file line number Diff line number Diff line change
Expand Up @@ -76,13 +76,11 @@ def __init__(self, outfilename, **kwargs):
self._runtyp = None

self.ts = self._Timestep(0) # need for properties initial calculations
self.periodic = False
self.delta = kwargs.pop("delta", 1.0) # there is no time, so let delta be 1

# update runtyp property
self.runtyp
if not self.runtyp in ['optimize', 'surface']:
raise AttributeError('Wrong RUNTYP= '+self.runtyp)
self.skip_timestep = 1

self.ts = self._Timestep(self.numatoms, **self._ts_kwargs)
# update numframes property
Expand Down Expand Up @@ -248,7 +246,6 @@ def open_trajectory(self):
# reset ts
ts = self.ts
ts.frame = -1
ts.time = 0
return self.outfile

def close(self):
Expand Down
9 changes: 5 additions & 4 deletions package/MDAnalysis/coordinates/INPCRD.py
Original file line number Diff line number Diff line change
Expand Up @@ -37,13 +37,14 @@ def _read_first_frame(self):
self.title = inf.readline().strip()
line = inf.readline().split()
self.numatoms = int(line[0])

self.ts = self._Timestep(self.numatoms, **self._ts_kwargs)
try:
time = float(line[1])
except IndexError:
time = 0.0

self.ts = self._Timestep(self.numatoms, **self._ts_kwargs)
self.ts.time = time
pass
else:
self.ts.time = time
self.ts.frame = 0

for p in xrange(self.numatoms // 2):
Expand Down
3 changes: 0 additions & 3 deletions package/MDAnalysis/coordinates/MOL2.py
Original file line number Diff line number Diff line change
Expand Up @@ -85,9 +85,6 @@ def __init__(self, filename, **kwargs):
self.substructure = {}
self.frames = blocks
self.numframes = len(blocks)
self.periodic = False
self.delta = 0
self.skip_timestep = 1

def parse_block(self, block):
sections = {}
Expand Down
12 changes: 0 additions & 12 deletions package/MDAnalysis/coordinates/PDB.py
Original file line number Diff line number Diff line change
Expand Up @@ -476,15 +476,6 @@ def __init__(self, filename, **kwargs):

frames = {}

ts_kwargs = {}
for att in ('dt', 'time_offset'):
try:
val = kwargs[att]
except KeyError:
pass
else:
ts_kwargs[att] = val

self.ts = self._Timestep(self._numatoms, **self._ts_kwargs)

pos = 0 # atom position for filling coordinates array
Expand Down Expand Up @@ -556,9 +547,6 @@ def __init__(self, filename, **kwargs):

self.frames = frames
self.numframes = len(frames) if frames else 1
self.periodic = False
self.delta = 0
self.skip_timestep = 1

def get_occupancy(self):
"""Return an array of occupancies in atom order."""
Expand Down
6 changes: 1 addition & 5 deletions package/MDAnalysis/coordinates/PDBQT.py
Original file line number Diff line number Diff line change
Expand Up @@ -176,11 +176,7 @@ def _c(start, stop, typeclass=float):
self.convert_pos_from_native(self.ts._pos) # in-place !
self.convert_pos_from_native(self.ts._unitcell[:3]) # in-place ! (only lengths)
self.numframes = 1
self.fixed = 0
self.skip = 1
self.periodic = False
self.delta = 0
self.skip_timestep = 1

# hack for PDBQTParser:
self._atoms = np.rec.fromrecords(atoms,
names="serial,name,resName,chainID,resSeq,occupancy,tempFactor,"
Expand Down
32 changes: 15 additions & 17 deletions package/MDAnalysis/coordinates/TRJ.py
Original file line number Diff line number Diff line change
Expand Up @@ -76,7 +76,7 @@
boxes.
* The trajectory does not contain time information so we simply set
the time step to 1 ps (or the user could provide it as kwarg *delta*)
the time step to 1 ps (or the user could provide it as kwarg *dt*)
* No direct access of frames is implemented, only iteration through
the trajectory.
Expand Down Expand Up @@ -177,7 +177,7 @@ class TRJReader(base.Reader):
`numatoms` is left at its default value of ``None``.
The length of a timestep is not stored in the trajectory itself but can
be set by passing the `delta` keyword argument to the constructor; it
be set by passing the `dt` keyword argument to the constructor; it
is assumed to be in ps. The default value is 1 ps.
Functionality is currently limited to simple iteration over the
Expand All @@ -187,6 +187,7 @@ class TRJReader(base.Reader):
.. versionchanged:: 0.11.0
Frames now 0-based instead of 1-based
kwarg 'delta' renamed to 'dt', for uniformity with other Readers
"""
format = 'TRJ'
units = {'time': 'ps', 'length': 'Angstrom'}
Expand All @@ -201,8 +202,6 @@ def __init__(self, filename, numatoms=None, **kwargs):
self._numframes = None

self.trjfile = None # have _read_next_timestep() open it properly!
self.skip_timestep = 1 # always 1 for trj at the moment
self.delta = kwargs.pop("delta", 1.0) # can set delta manually, default is 1ps
self.ts = self._Timestep(self.numatoms, **self._ts_kwargs)

# FORMAT(10F8.3) (X(i), Y(i), Z(i), i=1,NATOM)
Expand Down Expand Up @@ -351,7 +350,7 @@ def open_trajectory(self):
# reset ts
ts = self.ts
ts.frame = -1
ts.time = 0

return self.trjfile

def close(self):
Expand Down Expand Up @@ -403,6 +402,7 @@ class NCDFReader(base.Reader):
Added ability to read Forces
.. versionchanged:: 0.11.0
Frame labels now 0-based instead of 1-based
kwarg 'delta' renamed to 'dt', for uniformity with other Readers
"""

format = 'NCDF'
Expand Down Expand Up @@ -477,8 +477,6 @@ def __init__(self, filename, numatoms=None, **kwargs):
self.has_velocities = 'velocities' in self.trjfile.variables
self.has_forces = 'forces' in self.trjfile.variables

self.skip_timestep = 1 # always 1 for trj at the moment ? CHECK DOCS??
self.delta = kwargs.pop("delta", 1.0) # SHOULD GET FROM NCDF (as mean(diff(time)))??
self.periodic = 'cell_lengths' in self.trjfile.variables
self._current_frame = 0

Expand Down Expand Up @@ -549,7 +547,7 @@ def Writer(self, filename, **kwargs):
:Keywords:
*numatoms*
number of atoms
*delta*
*dt*
length of one timestep in picoseconds
*remarks*
string that is stored in the title field
Expand All @@ -558,7 +556,7 @@ def Writer(self, filename, **kwargs):
"""
numatoms = kwargs.pop('numatoms', self.numatoms)
kwargs.setdefault('remarks', self.remarks)
kwargs.setdefault('delta', self.delta)
kwargs.setdefault('dt', self.dt)
return NCDFWriter(filename, numatoms, **kwargs)


Expand All @@ -583,14 +581,16 @@ class NCDFWriter(base.Writer):
.. versionchanged:: 0.10.0
Added ability to write velocities and forces
.. versionchanged:: 0.11.0
kwarg 'delta' renamed to 'dt', for uniformity with other Readers
"""

format = 'NCDF'
version = "1.0"
units = {'time': 'ps', 'length': 'Angstrom', 'velocity': 'Angstrom/ps',
'force': 'kcal/(mol*Angstrom)'}

def __init__(self, filename, numatoms, start=0, step=1, delta=1.0, remarks=None,
def __init__(self, filename, numatoms, start=0, step=1, dt=1.0, remarks=None,
convert_units=None, zlib=False, cmplevel=1, **kwargs):
"""Create a new NCDFWriter
Expand All @@ -605,7 +605,7 @@ def __init__(self, filename, numatoms, start=0, step=1, delta=1.0, remarks=None,
starting timestep
*step*
skip between subsequent timesteps
*delta*
*dt*
timestep
*convert_units*
``True``: units are converted to the AMBER base format; ``None`` selects
Expand All @@ -630,7 +630,7 @@ def __init__(self, filename, numatoms, start=0, step=1, delta=1.0, remarks=None,

self.start = start # do we use those?
self.step = step # do we use those?
self.delta = delta
self.dt = dt
self.remarks = remarks or "AMBER NetCDF format (MDAnalysis.coordinates.trj.NCDFWriter)"

self.zlib = zlib
Expand Down Expand Up @@ -771,14 +771,12 @@ def _write_next_timestep(self, ts):
try:
time = self.convert_time_to_native(ts.time, inplace=False)
except AttributeError:
time = ts.frame * self.convert_time_to_native(self.delta, inplace=False)
time = ts.frame * self.convert_time_to_native(self.dt, inplace=False)
unitcell = self.convert_dimensions_to_unitcell(ts)
else:
pos = ts._pos
try:
time = ts.time
except AttributeError:
time = ts.frame * self.delta
time = ts.time

unitcell = ts.dimensions

# write step
Expand Down

0 comments on commit 0e36461

Please sign in to comment.