pytraj.TrajectoryIterator

out-of-core TrajectoryIterator

class pytraj.trajectory_iterator.TrajectoryIterator(filename=None, top=None, *args, **kwd)

Bases: pytraj.c_traj.c_trajectory.TrajectoryCpptraj

out-of-core trajectory holder.

Notes

It’s a bit tricky to pickle this class. As default, new TrajectoryIterator will use original trj filename and top filename. If set _pickle_topology to True, its Topology will be pickled (slow but more accurate if you change the topology in the fly)

Examples

>>> import pytraj as pt
>>> from pytraj.testing import get_fn
>>> traj_name, top_name = get_fn('tz2')
>>> traj = pt.TrajectoryIterator(traj_name, top_name)
>>> # user should always use :method:`pytraj.iterload` to load TrajectoryIterator
>>> traj = pt.iterload(['remd.x.000', 'remd.x.001'], 'input.parm7') 

Attributes

filelist return a list of input filenames. Order does matter
filename return 1st filename in filelist. For testing only
metadata return a dict of general information
n_atoms used for frame_iter
n_frames
shape (n_frames, n_atoms, 3)
temperatures
top
topology traditional name for Topology file
unitcells return unitcells (Box) with shape=(n_frames, 6)
xyz return 3D array of coordinates

Methods

__call__ Call self as a function.
at(index)
copy() return a deep copy. Use this method with care since the copied traj just reuse
iterchunk([chunksize, start, stop, ...]) iterate trajectory by chunk
iterframe([start, stop, step, mask, ...]) iterate trajectory with given frame_indices or given (start, stop, step)
save(self[, filename, overwrite]) convenient method to save Trajectory
at(index)
copy()

return a deep copy. Use this method with care since the copied traj just reuse the filenames

filelist

return a list of input filenames. Order does matter

filename

return 1st filename in filelist. For testing only

iterchunk(chunksize=2, start=0, stop=-1, autoimage=False, rmsfit=None)

iterate trajectory by chunk

Parameters:

chunk : int, default=2

size of each chunk. Notes: final chunk’s size might be changed

start : int, default=0 (first frame)

stop : int, default=-1 (last frame)

autoimage : bool, default=False

rmsfit : None | tuple/list of (reference frame, mask)

Notes

if using ‘autoimage` with reference frame for rms-fit, make sure to autoimage ref first

Examples

>>> import pytraj as pt
>>> traj = pt.load_sample_data('tz2')
>>> ref = traj[3]
>>> for chunk in traj.iterchunk(3, autoimage=True, rmsfit=(ref, '@CA')): pass
iterframe(start=0, stop=None, step=1, mask=None, autoimage=False, rmsfit=None, copy=False, frame_indices=None)

iterate trajectory with given frame_indices or given (start, stop, step)

Parameters:

start : int, default 0

stop : {None, int}, default None

if None, iterate to final frame

step : int, default 1

mask : {None, str}, default None

if None, use all atoms. If not None, use given mask

autoimage : bool, default False

if True, perform autoimage for each frame

rmsfit : {None, int, tuple}, default None

if not None, perform superpose each Frame to to reference.

frame_indices : {None, array-like}

if not None, iterate trajectory for given indices. If frame_indices is given, (start, stop, step) will be ignored.

Examples

>>> import pytraj as pt
>>> traj = pt.load_sample_data('tz2')
>>> for frame in traj.iterframe(0, 8, 2): pass
>>> for frame in traj.iterframe(0, 8, 2, autoimage=True): pass
>>> # use negative index
>>> traj.n_frames
10
>>> fi = traj.iterframe(0, -1, 2, autoimage=True)
>>> fi.n_frames
5
>>> # mask is atom indices
>>> fi = traj.iterframe(0, -1, 2, mask=range(100), autoimage=True)
>>> fi.n_atoms
100
metadata

return a dict of general information

Examples

>>> traj.metadata
{'box_type': 'ortho',
 'has_box': True,
 'has_force': False,
 'has_replcica_dims': False,
 'has_temperature': False,
 'has_time': True,
 'has_velocity': False,
 'n_atoms': 5293,
 'n_frames': 10}
n_atoms

used for frame_iter

n_frames
save(self, filename='', overwrite=True, **kwd)

convenient method to save Trajectory

Examples

>>> # save to DCD file
>>> traj.save('traj.dcd')
>>> # save to AMBER netcdf file
>>> traj.save('traj.nc')
shape

(n_frames, n_atoms, 3)

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> traj.shape
(10, 5293, 3)
temperatures
top
topology

traditional name for Topology file

Examples

>>> import pytraj as pt
>>> from pytraj.testing import get_fn
>>> fname, tname = get_fn('ala3')
>>> traj = pt.iterload(fname, tname)
>>> traj.topology
<Topology: 34 atoms, 3 residues, 1 mols, non-PBC>
>>> new_traj = pt.TrajectoryIterator()
>>> new_traj.topology = traj.topology
>>> new_traj._load(traj.filename)
unitcells

return unitcells (Box) with shape=(n_frames, 6)

xyz

return 3D array of coordinates