NetCDFFile¶
- class mdhelper.openmm.file.NetCDFFile(file: str | Dataset, mode: str, restart: bool = False, **kwargs)[source]¶
Bases:
object
Interface for reading and writing AMBER NetCDF trajectory and restart files.
- Parameters:
- filestr or netcdf4.Dataset
NetCDF file. If file is a filename and does not have the
.nc
or.ncdf
extension,.nc
will automatically be appended.- modestr
NetCDF file access mode.
- restartbool, default:
False
Specifies whether the NetCDF file is a trajectory or restart file.
- **kwargs
Keyword arguments to be passed to
netCDF4.Dataset
.
Methods
Get the simulation box dimensions.
Get the forces acting on the atoms.
Get the number of atoms.
Get the number of frames.
Get the atom positions.
Get simulation times.
Get atom velocities.
Write a single simulation state to a restart NetCDF file.
Initialize a NetCDF file according to AMBER NetCDF Trajectory/Restart Convention Version 1.0, Revision C.
Write simulation state(s) to a NetCDF file.
- get_dimensions(frames: int | list[int] | slice = None, units: bool = True) tuple[ndarray[float], ndarray[float]] | tuple[Quantity, Quantity] [source]¶
Get the simulation box dimensions.
- Parameters:
- framesint, list, or slice, optional
Frame indices. If
None
, the dimensions across all frames are returned.- unitsbool, default:
True
Determines whether the dimensions are returned with units.
- Returns:
- cell_lengthsnumpy.ndarray, optional
Simulation box dimensions.
Reference unit: \(\mathrm{Å}\).
- cell_anglesnumpy.ndarray, optional
Angles that define the shape of the simulation box.
Reference unit: \(^\circ\).
- get_forces(frames: int | list[int] | slice = None, units: bool = True) ndarray[float] | Quantity [source]¶
Get the forces acting on the atoms.
- Parameters:
- framesint, list, or slice, optional
Frame indices. If
None
, the forces across all frames are returned.- unitsbool, default:
True
Determines whether the forces are returned with units.
- Returns:
- forcesnumpy.ndarray or openmm.unit.Quantity
Forces acting on the atoms. If the NetCDF file does not contain this information,
None
is returned.Reference unit: \(\mathrm{Å/ps}\).
- get_positions(frames: int | list[int] | slice = None, units: bool = True) ndarray[float] | Quantity [source]¶
Get the atom positions.
- Parameters:
- framesint, list, or slice, optional
Frame indices. If
None
, the positions across all frames are returned.- unitsbool, default:
True
Determines whether the positions are returned with units.
- Returns:
- positionsnumpy.ndarray or openmm.unit.Quantity
Atom positions.
Reference unit: \(\mathrm{Å}\).
- get_times(frames: int | list[int] | slice = None, units: bool = True) ndarray[float] | Quantity [source]¶
Get simulation times.
- Parameters:
- framesint, list, or slice, optional
Frame indices. If
None
, the times across all frames are returned.- unitsbool, default:
True
Determines whether the times are returned with units.
- Returns:
- timesnumpy.ndarray or openmm.unit.Quantity
Simulation times.
Reference unit: \(\mathrm{ps}\).
- get_velocities(frames: int | list[int] | slice = None, units: bool = True) ndarray[float] | Quantity [source]¶
Get atom velocities.
- Parameters:
- framesint, list, or slice, optional
Frame indices. If
None
, the velocities across all frames are returned.- unitsbool, default:
True
Determines whether the velocities are returned with units.
- Returns:
- velocitiesnumpy.ndarray or openmm.unit.Quantity
Atom velocities. If the NetCDF file does not contain this information,
None
is returned.Reference unit: \(\mathrm{Å/ps}\).
- write_file(state: State) NetCDFFile [source]¶
Write a single simulation state to a restart NetCDF file.
Note
This function can be used as either a static or instance method.
- Parameters:
- selfstr, netcdf4.Dataset, or mdhelper.openmm.file.NetCDFFile
If
write_file()
is called as a static method, you must provide a filename or a NetCDF file object. Otherwise, the NetCDF file embedded in the current instance is used.- stateopenmm.State
OpenMM simulation state from which to retrieve cell dimensions and atom positions, velocities, and forces.
- Returns:
- netcdf_filemdhelper.openmm.file.NetCDFFile
NetCDF file object. Only returned when this function is used as a static method.
- write_header(N: int, cell: bool, velocities: bool, forces: bool, restart: bool = False, *, remd: str = None, temp0: float = None, remd_dimtype: ndarray[int] = None, remd_indices: ndarray[int] = None, remd_repidx: int = -1, remd_crdidx: int = -1, remd_values: ndarray[float] = None) NetCDFFile [source]¶
Initialize a NetCDF file according to AMBER NetCDF Trajectory/Restart Convention Version 1.0, Revision C.
Note
This function can be used as either a static or instance method.
- Parameters:
- selfstr, netcdf4.Dataset, or mdhelper.openmm.file.NetCDFFile
If
write_header()
is called as a static method, you must provide a filename or a NetCDF file object. Otherwise, the NetCDF file embedded in the current instance is used.- Nint
Number of atoms.
- cellbool
Specifies whether simulation box length and angle information is available.
- velocitiesbool
Specifies whether atom velocities should be written.
- forcesbool
Specifies whether forces exerted on atoms should be written.
- restartbool, default:
False
Specifies whether the NetCDF file is a trajectory or restart file.
- remdstr, keyword-only, optional
Specifies whether information about a replica exchange molecular dynamics (REMD) simulation is written.
Valid values:
"temp"
for regular REMD."multi"
for multi-dimensional REMD.
- temp0float, keyword-only, optional
Temperature that the thermostat is set to maintain for a REMD restart file only.
Reference unit: \(\mathrm{K}\).
- remd_dimtypearray-like, keyword-only, optional
Array specifying the exchange type(s) for the REMD dimension(s). Required for a multi-dimensional REMD restart file.
- remd_indicesarray-like, keyword-only, optional
Array specifying the position in all dimensions that each frame is in. Required for a multi-dimensional REMD restart file.
- remd_repidxint, keyword-only, optional
Overall index of the frame in replica space.
- remd_crdidxint, keyword-only, optional
Overall index of the frame in coordinate space.
- remd_valuesarray-like, keyword-only, optional
Replica value the specified replica dimension has for that given frame. Required for a multi-dimensional REMD restart file.
- Returns:
- netcdf_filemdhelper.openmm.file.NetCDFFile
NetCDF file object. Only returned when this function is used as a static method.
- write_model(time: float | ndarray[float], coordinates: ndarray[float], velocities: ndarray[float] = None, forces: ndarray[float] = None, cell_lengths: ndarray[float] = None, cell_angles: ndarray[float] = None, *, restart: bool = False) NetCDFFile [source]¶
Write simulation state(s) to a NetCDF file.
Note
This function can be used as either a static or instance method.
- Parameters:
- selfstr, netcdf4.Dataset, or mdhelper.openmm.file.NetCDFFile
If
write_model()
is called as a static method, you must provide a filename or a NetCDF file object. Otherwise, the NetCDF file embedded in the current instance is used.- timefloat or numpy.ndarray
Time(s). The dimensionality determines whether a single or multiple frames are written.
Reference unit: \(\mathrm{ps}\).
- coordinatesnumpy.ndarray
Coordinates of \(N\) atoms over \(N_t\) frames. The dimensionality depends on whether a single or multiple frames are to be written and must be compatible with that for time.
Shape: \((N,\,3)\) or \((N_t,\,N,\,3)\).
Reference unit: \(\mathrm{Å}\).
- velocitiesnumpy.ndarray, optional
Velocities of \(N\) atoms over \(N_t\) frames. The dimensionality depends on whether a single or multiple frames are to be written and must be compatible with that for time.
Shape: \((N,\,3)\) or \((N_t,\,N,\,3)\).
Reference unit: \(\mathrm{Å/ps}\).
- forcesnumpy.ndarray, optional
Forces exerted on \(N\) atoms over \(N_t\) frames. The dimensionality depends on whether a single or multiple frames are to be written and must be compatible with that for time.
Shape: \((N,\,3)\) or \((N_t,\,N,\,3)\).
Reference unit: \(\mathrm{Å/ps}\).
- cell_lengthsnumpy.ndarray, optional
Simulation box dimensions.
Shape: \((3,)\).
Reference unit: \(\mathrm{Å}\).
- cell_anglesnumpy.ndarray, optional
Angles that define the shape of the simulation box.
Shape: \((3,)\).
Reference unit: \(^\circ\).
- restartbool, keyword-only, default:
False
Prevents the frame index from being incremented if writing a NetCDF restart file.
- Returns:
- netcdf_filemdhelper.openmm.file.NetCDFFile
NetCDF file object. Only returned when this function is used as a static method.