Most commmon usage: load a single file with all frame
# read amber format
In [1]: import pytraj as pt
In [2]: traj = pt.iterload('tz2.ortho.nc', 'tz2.ortho.parm7')
In [3]: traj
Out[3]: 
pytraj.TrajectoryIterator, 10 frames: 
Size: 0.001183 (GB)
<Topology: 5293 atoms, 1704 residues, 1692 mols, PBC with box type = ortho>
Most commmon usage: load a single file with frame stride
# load only frame 0, 2, 4, 6 (python use 0-based index and skip final index)
traj = pt.iterload('tz2.nc', 'tz2.parm7', frame_slice=[(0, 8, 2),])
Load many files at once
In [1]: ls remd.x.*
remd.x.000  remd.x.003  remd.x.006  remd.x.009  remd.x.012  remd.x.015
remd.x.001  remd.x.004  remd.x.007  remd.x.010  remd.x.013  remd.x.016
remd.x.002  remd.x.005  remd.x.008  remd.x.011  remd.x.014  remd.x.017
In [2]: traj = pt.iterload('remd.x.*', 'myparm.parm7')
In [3]: traj
Out[3]:
<pytraj.TrajectoryIterator, 40000 frames, include:
<Topology: 17443 atoms, 5666 residues, 5634 mols, PBC with box type = truncoct>>
In [4]: traj._estimated_MB
Out[4]: 15969.54345703125
Load many files with frame stride
Example: load 1000 frames from two trajectories (500 each), skip every two frames.
In [1]: traj = pt.iterload(['remd.x.000', 'remd.x.001'], 'myparm.parm7', frame_slice=[(0, 1000, 2), (0, 1000, 2)])
In [2]: traj
Out[2]:
<pytraj.TrajectoryIterator, 1000 frames, include:
<Topology: 17443 atoms, 5666 residues, 5634 mols, PBC with box type = truncoct>>
which is similiar to cpptraj input:
parm myparm.parm7
trajin remd.x.000 1 1000 2
trajin remd.x.001 1 1000 2
Note
cpptraj uses 1-based index for input while pytraj used 0-based index.
# write to CHARMM dcd format
In [4]: pt.write_traj('test.dcd', traj, overwrite=True)
# write with given frames
In [5]: pt.write_traj('test2.dcd', traj, frame_indices=[0, 3, 7, 9], overwrite=True)
pytraj.write_traj(filename, traj, format='infer', top=None, frame_indices=None, overwrite=False, force=False, velocity=False, options='')¶write Trajectory-like or iterable object to trajectory file
| Parameters: | filename : str traj : Trajectory-like or iterator that produces Frame or 3D ndarray with shape=(n_frames, n_atoms, 3) format : str, default ‘infer’ 
 top : Topology, optional, default: None frame_indices: array-like or iterator that produces integer, default: None 
 overwrite: bool, default: False velocity : bool, default False 
 force : bool, default False 
 options : str, additional cpptraj keywords | 
|---|
Notes
| Format | Extension | 
|---|---|
| Amber Trajectory | .crd | 
| Amber NetCDF | .nc | 
| Amber Restart | .rst7 | 
| Amber NetCDF | .ncrst | 
| Charmm DCD | .dcd | 
| PDB | .pdb | 
| Mol2 | .mol2 | 
| Scripps | .binpos | 
| Gromacs | .trr | 
| SQM Input | .sqm | 
‘options’ for writing to pdb format (cptraj manual):
dumpq:       Write atom charge/GB radius in occupancy/B-factor columns (PQR format)."
parse:       Write atom charge/PARSE radius in occupancy/B-factor columns (PQR format)."
vdw:         Write atom charge/VDW radius in occupancy/B-factor columns (PQR format)."
pdbres:      Use PDB V3 residue names."
pdbatom:     Use PDB V3 atom names."
pdbv3:       Use PDB V3 residue/atom names."
teradvance:  Increment record (atom) # for TER records (default no)."
terbyres:    Print TER cards based on residue sequence instead of molecules."
model:       Write to single file separated by MODEL records."
multi:       Write each frame to separate files."
chainid <c>: Write character 'c' in chain ID column."
sg <group>:  Space group for CRYST1 record, only if box coordinates written."
include_ep:  Include extra points."
conect:      Write CONECT records using bond information.");
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> pt.write_traj("output/t.nc", traj, overwrite=True) # write to amber netcdf file
>>> # write to multi pdb files (t.pdb.1, t.pdb.2, ...)
>>> pt.write_traj("output/t.pdb", traj, overwrite=True, options='multi')
>>> # write all frames to single pdb file and each frame is seperated by "MODEL" word
>>> pt.write_traj("output/t.pdb", traj, overwrite=True, options='model')
>>> # write to DCD file
>>> pt.write_traj("output/test.dcd", traj, overwrite=True)
>>> # write to netcdf file from 3D numpy array, need to provide Topology
>>> xyz = traj.xyz
>>> top = traj.top
>>> pt.write_traj("output/test_xyz.nc", xyz, top=traj.top, overwrite=True)
>>> pt.write_traj("output/test_xyz.nc", xyz, top=traj.top, overwrite=True)
>>> # you can make a fake Topology to write xyz coordinates too
>>> n_atoms = xyz.shape[1]
>>> top2 = pt.tools.make_fake_topology(n_atoms)
>>> pt.write_traj("output/test_xyz_fake_top.nc", xyz, top=top2, overwrite=True)
‘options’ for writing to amber netcdf format (cptraj manual):
remdtraj: Write temperature to trajectory (makes REMD trajectory)."
velocity: Write velocities to trajectory."
force: Write forces to trajectory.");
‘options’ for writing to amber netcdf restart format(cptraj manual):
novelocity: Do not write velocities to restart file."
notime:     Do not write time to restart file."
remdtraj:   Write temperature to restart file."
time0:      Time for first frame (default 1.0)."
dt:         Time step for subsequent frames, t=(time0+frame)*dt; (default 1.0)");
keepext     Keep filename extension; write '<name>.<num>.<ext>' instead (example: myfile.1.rst7)
‘options’ for writing to mol2 format (cptraj manual):
single   : Write to a single file."
multi    : Write each frame to a separate file."
sybyltype: Convert Amber atom types (if present) to SYBYL types.");
‘options’ for other formats:
please check http://ambermd.org/doc12/Amber15.pdf