6.20. TRZ trajectory I/O  — MDAnalysis.coordinates.TRZ
Classes to read IBIsCO / YASP TRZ binary trajectories, including
coordinates, velocities and more (see attributes of the Timestep).
Data are read and written in binary representation but because this depends on the machine hardware architecture, MDAnalysis always reads and writes TRZ trajectories in little-endian byte order.
6.20.1. Classes
- class MDAnalysis.coordinates.TRZ.TRZReader(trzfilename, n_atoms=None, **kwargs)[source]
- Reads an IBIsCO or YASP trajectory file - Note - Binary TRZ trajectories are always assumed to be written in little-endian byte order and are read as such. - Changed in version 0.11.0: Frames now 0-based instead of 1-based. Extra data (Temperature, Energies, Pressures, etc) now read into ts.data dictionary. Now passes a weakref of self to ts (ts._reader). - Changed in version 1.0.1: Now checks for the correct n_atoms on reading and can raise - ValueError.- Changed in version 2.1.0: TRZReader now returns a default - dtof 1.0 when it cannot be obtained from the difference between two frames.- Creates a TRZ Reader - Parameters
- Raises
- ValueError – If n_atoms or the number of atoms in the topology file do not match the number of atoms in the trajectory. 
 - Writer(filename, n_atoms=None)[source]
- A trajectory writer with the same properties as this trajectory. 
 - property delta
- MD integration timestep 
 - property n_atoms
- Number of atoms in a frame 
 - property n_frames
- Total number of frames in a trajectory 
 - property skip_timestep
- Timesteps between trajectory frames 
 - units = {'length': 'nm', 'time': 'ps', 'velocity': 'nm/ps'}
- dict with units of of time and length (and velocity, force, … for formats that support it) 
 
- class MDAnalysis.coordinates.TRZ.TRZWriter(filename, n_atoms, title='TRZ', convert_units=True)[source]
- Writes a TRZ format trajectory. - Note - Binary TRZ trajectories are always written in little-endian byte order. - Create a TRZWriter - Parameters
- filename (str) – name of output file 
- n_atoms (int) – number of atoms in trajectory 
- title (str (optional)) – title of the trajectory; the title must be 80 characters or shorter, a longer title raises a ValueError exception. 
- convert_units (bool (optional)) – units are converted to the MDAnalysis base format; [ - True]
 
 - units = {'length': 'nm', 'time': 'ps', 'velocity': 'nm/ps'}
- dict with units of of time and length (and velocity, force, … for formats that support it)