mdtraj.formats.NetCDFTrajectoryFile¶
-
class
mdtraj.formats.
NetCDFTrajectoryFile
(filename, mode='r', force_overwrite=True)¶ Interface for reading and writing to AMBER NetCDF files. This is a file-like object, that supports both reading or writing depending on the mode flag. It implements the context manager protocol, so you can also use it with the python ‘with’ statement.
- Parameters
- filenamestr
The name of the file to open
- mode{‘r’, ‘w’}, default=’r’
The mode in which to open the file. Valid options are ‘r’ and ‘w’ for ‘read’ and ‘write’, respectively.
- force_overwritebool, default=False
In write mode, if a file named filename already exists, clobber it and overwrite it.
- Attributes
- n_atoms
- n_frames
Methods
close
(self)Close the NetCDF file handle
flush
(self)Write all buffered data in the to the disk file.
read
(self[, n_frames, stride, atom_indices])Read data from a molecular dynamics trajectory in the AMBER NetCDF format.
read_as_traj
(self, topology[, n_frames, …])Read a trajectory from a NetCDF file
seek
(self, offset[, whence])Move to a new file position
tell
(self)Current file position
write
(self, coordinates[, time, …])Write one or more frames of a molecular dynamics trajectory to disk in the AMBER NetCDF format.
-
__init__
(self, filename, mode='r', force_overwrite=True)¶ Initialize self. See help(type(self)) for accurate signature.
Methods
__init__
(self, filename[, mode, force_overwrite])Initialize self.
close
(self)Close the NetCDF file handle
flush
(self)Write all buffered data in the to the disk file.
read
(self[, n_frames, stride, atom_indices])Read data from a molecular dynamics trajectory in the AMBER NetCDF format.
read_as_traj
(self, topology[, n_frames, …])Read a trajectory from a NetCDF file
seek
(self, offset[, whence])Move to a new file position
tell
(self)Current file position
write
(self, coordinates[, time, …])Write one or more frames of a molecular dynamics trajectory to disk in the AMBER NetCDF format.
Attributes
distance_unit
n_atoms
n_frames
-
close
(self)¶ Close the NetCDF file handle
-
flush
(self)¶ Write all buffered data in the to the disk file.
-
read
(self, n_frames=None, stride=None, atom_indices=None)¶ Read data from a molecular dynamics trajectory in the AMBER NetCDF format.
- Parameters
- n_framesint, optional
If n_frames is not None, the next n_frames of data from the file will be read. Otherwise, all of the frames in the file will be read.
- strideint, optional
If stride is not None, read only every stride-th frame from disk.
- atom_indicesnp.ndarray, dtype=int, optional
The specific indices of the atoms you’d like to retrieve. If not supplied, all of the atoms will be retrieved.
- Returns
- coordinatesnp.ndarray, shape=(n_frames, n_atoms, 3)
The cartesian coordinates of the atoms, in units of angstroms.
- timenp.ndarray, None
The time corresponding to each frame, in units of picoseconds, or None if no time information is present in the trajectory.
- cell_lengthsnp.ndarray, None
The lengths (a,b,c) of the unit cell for each frame, or None if the information is not present in the file.
- cell_anglesnp.ndarray, None
The angles (lpha, eta, gamma) defining the unit cell for each frame, or None if the information is not present in the file.
-
read_as_traj
(self, topology, n_frames=None, stride=None, atom_indices=None)¶ Read a trajectory from a NetCDF file
- Parameters
- topologyTopology
The system topology
- n_framesint, optional
If positive, then read only the next n_frames frames. Otherwise read all of the frames in the file.
- stridenp.ndarray, optional
Read only every stride-th frame.
- atom_indicesarray_like, optional
If not none, then read only a subset of the atoms coordinates from the file. This may be slightly slower than the standard read because it required an extra copy, but will save memory.
- Returns
- trajectoryTrajectory
A trajectory object containing the loaded portion of the file.
-
seek
(self, offset, whence=0)¶ Move to a new file position
- Parameters
- offsetint
A number of frames.
- whence{0, 1, 2}
0: offset from start of file, offset should be >=0. 1: move relative to the current position, positive or negative 2: move relative to the end of file, offset should be <= 0. Seeking beyond the end of a file is not supported
-
tell
(self)¶ Current file position
- Returns
- offsetint
The current frame in the file.
-
write
(self, coordinates, time=None, cell_lengths=None, cell_angles=None)¶ Write one or more frames of a molecular dynamics trajectory to disk in the AMBER NetCDF format.
- Parameters
- coordinatesnp.ndarray, dtype=np.float32, shape=(n_frames, n_atoms, 3)
The cartesian coordinates of each atom, in units of angstroms.
- timenp.ndarray, dtype=np.float32, shape=(n_frames), optional
The time index corresponding to each frame, in units of picoseconds.
- cell_lengthsnp.ndarray, dtype=np.double, shape=(n_frames, 3)
The lengths (a,b,c) of the unit cell for each frame.
- cell_anglesnp.ndarray, dtype=np.double, shape=(n_frames, 3)
The angles (lpha, eta, gamma) defining the unit cell for each frame.
Notes
If the input arrays are of dimension deficient by one, for example if the coordinates array is two dimensional, the time is a single scalar or cell_lengths and cell_angles are a 1d array of length three, that is okay. You’ll simply be saving a single frame.