7.11. LAMMPS DCD trajectory and DATA I/O — MDAnalysis.coordinates.LAMMPS
Classes to read and write LAMMPS DCD binary trajectories, LAMMPS DATA files and LAMMPS dump files. Trajectories can be read regardless of system-endianness as this is auto-detected.
LAMMPS can write DCD trajectories but unlike a CHARMM trajectory (which is often called a DCD even though CHARMM itself calls them “trj”) the time unit is not fixed to be the AKMA time unit (20 AKMA is 0.978 picoseconds or 1 AKMA = 4.888821e-14 s) but can depend on settings in LAMMPS. The most common case for biomolecular simulations appears to be that the time step is recorded in femtoseconds (command units real in the input file) and lengths in ångströms. Other cases are unit-less Lennard-Jones time units.
This presents a problem for MDAnalysis because it cannot autodetect
the unit from the file. By default we are assuming that the unit for
length is the ångström and for the time is the femtosecond. If this is
not true then the user should supply the appropriate units in the
keywords timeunit and/or lengthunit to DCDWriter
and
Universe
(which then calls
DCDReader
).
7.11.1. Data file formats
By default either the atomic or full atom styles are expected, however this can be customised, see Atom styles.
7.11.2. Dump files
The DumpReader expects ascii dump files written with the default LAMMPS dump format of ‘atom’
7.11.3. Example: Loading a LAMMPS simulation
To load a LAMMPS simulation from a LAMMPS data file (using the
DATAParser
) together with a
LAMMPS DCD with “real” provide the keyword format=”LAMMPS”:
>>> import MDAnalysis
>>> from MDAnalysis.tests.datafiles import LAMMPSdata2, LAMMPSdcd2
>>> u = MDAnalysis.Universe(LAMMPSdata2, LAMMPSdcd2, format="LAMMPS")
If the trajectory uses units nano then use
>>> import MDAnalysis
>>> from MDAnalysis.tests.datafiles import LAMMPSdata2, LAMMPSdcd2
>>> u = MDAnalysis.Universe(LAMMPSdata2, LAMMPSdcd2, format="LAMMPS",
... lengthunit="nm", timeunit="ns")
To scan through a trajectory to find a desirable frame and write to a LAMMPS data file,
>>> import MDAnalysis
>>> from MDAnalysis.tests.datafiles import LAMMPSdata2, LAMMPSdcd2
>>> u = MDAnalysis.Universe(LAMMPSdata2, LAMMPSdcd2, format="LAMMPS",
... lengthunit="nm", timeunit="ns")
>>> take_this_frame = False
>>> for ts in u.trajectory:
... # analyze frame
... if ts.frame == 4:
... take_this_frame = True
... if take_this_frame == True:
... with MDAnalysis.Writer('frame.data') as W:
... W.write(u.atoms)
... break
Note
Lennard-Jones units are not implemented. See MDAnalysis.units
for other recognized values and the documentation for the LAMMPS
units command.
7.11.4. Classes
- class MDAnalysis.coordinates.LAMMPS.DCDReader(dcdfilename, **kwargs)[source]
Read a LAMMPS DCD trajectory.
The units can be set from the constructor with the keyword arguments timeunit and lengthunit. The defaults are “fs” and “Angstrom”, corresponding to LAMMPS units style “real”. See
MDAnalysis.units
for other recognized values.- Parameters:
Changed in version 0.17.0: Changed to use libdcd.pyx library and removed the correl function
- OtherWriter(filename, **kwargs)
Returns a writer appropriate for filename.
Sets the default keywords start, step and dt (if available). n_atoms is always set from
Reader.n_atoms
.See also
Reader.Writer()
- add_auxiliary(aux_spec: str | Dict[str, str] | None = None, auxdata: str | AuxReader | None = None, format: str | None = None, **kwargs) None
Add auxiliary data to be read alongside trajectory.
Auxiliary data may be any data timeseries from the trajectory additional to that read in by the trajectory reader. auxdata can be an
AuxReader
instance, or the data itself as e.g. a filename; in the latter case an appropriateAuxReader
is guessed from the data/file format. An appropriate format may also be directly provided as a key word argument.On adding, the AuxReader is initially matched to the current timestep of the trajectory, and will be updated when the trajectory timestep changes (through a call to
next()
or jumping timesteps withtrajectory[i]
).The representative value(s) of the auxiliary data for each timestep (as calculated by the
AuxReader
) are stored in the current timestep in thets.aux
namespace under aux_spec; e.g. to add additional pull force data stored in pull-force.xvg:u = MDAnalysis.Universe(PDB, XTC) u.trajectory.add_auxiliary('pull', 'pull-force.xvg')
The representative value for the current timestep may then be accessed as
u.trajectory.ts.aux.pull
oru.trajectory.ts.aux['pull']
.The following applies to energy readers like the
EDRReader
.All data that is present in the (energy) file can be added by omitting aux_spec like so:
u.trajectory.add_auxiliary(auxdata="ener.edr")
aux_spec is expected to be a dictionary that maps the desired attribute name in the
ts.aux
namespace to the precise data to be added as identified by adata_selector
:term_dict = {"temp": "Temperature", "epot": "Potential"} u.trajectory.add_auxiliary(term_dict, "ener.edr")
Adding this data can be useful, for example, to filter trajectory frames based on non-coordinate data like the potential energy of each time step. Trajectory slicing allows working on a subset of frames:
selected_frames = np.array([ts.frame for ts in u.trajectory if ts.aux.epot < some_threshold]) subset = u.trajectory[selected_frames]
See also
Note
Auxiliary data is assumed to be time-ordered, with no duplicates. See the Auxiliary API.
- add_transformations(*transformations)
Add all transformations to be applied to the trajectory.
This function take as list of transformations as an argument. These transformations are functions that will be called by the Reader and given a
Timestep
object as argument, which will be transformed and returned to the Reader. The transformations can be part of thetransformations
module, or created by the user, and are stored as a list transformations. This list can only be modified once, and further calls of this function will raise an exception.u = MDAnalysis.Universe(topology, coordinates) workflow = [some_transform, another_transform, this_transform] u.trajectory.add_transformations(*workflow)
The transformations are applied in the order given in the list transformations, i.e., the first transformation is the first or innermost one to be applied to the
Timestep
. The example above would be equivalent tofor ts in u.trajectory: ts = this_transform(another_transform(some_transform(ts)))
- Parameters:
transform_list (list) – list of all the transformations that will be applied to the coordinates in the order given in the list
See also
- property aux_list
Lists the names of added auxiliary data.
- check_slice_indices(start, stop, step)
Check frame indices are valid and clip to fit trajectory.
The usage follows standard Python conventions for
range()
but see the warning below.- Parameters:
start (int or None) – Starting frame index (inclusive).
None
corresponds to the default of 0, i.e., the initial frame.stop (int or None) – Last frame index (exclusive).
None
corresponds to the default of n_frames, i.e., it includes the last frame of the trajectory.step (int or None) – step size of the slice,
None
corresponds to the default of 1, i.e, include every frame in the range start, stop.
- Returns:
start, stop, step – Integers representing the slice
- Return type:
Warning
The returned values start, stop and step give the expected result when passed in
range()
but gives unexpected behavior when passed in aslice
whenstop=None
andstep=-1
This can be a problem for downstream processing of the output from this method. For example, slicing of trajectories is implemented by passing the values returned by
check_slice_indices()
torange()
range(start, stop, step)
and using them as the indices to randomly seek to. On the other hand, in
MDAnalysis.analysis.base.AnalysisBase
the values returned bycheck_slice_indices()
are used to splice the trajectory by creating aslice
instanceslice(start, stop, step)
This creates a discrepancy because these two lines are not equivalent:
range(10, -1, -1) # [10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0] range(10)[slice(10, -1, -1)] # []
- convert_forces_from_native(force, inplace=True)
Conversion of forces array force from native to base units
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_forces_to_native(force, inplace=True)
Conversion of force array force from base to native units.
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_pos_from_native(x, inplace=True)
Conversion of coordinate array x from native units to base units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_pos_to_native(x, inplace=True)
Conversion of coordinate array x from base units to native units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_from_native(t, inplace=True)
Convert time t from native units to base units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.) In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_to_native(t, inplace=True)
Convert time t from base units to native units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool, optional) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_velocities_from_native(v, inplace=True)
Conversion of velocities array v from native to base units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- convert_velocities_to_native(v, inplace=True)
Conversion of coordinate array v from base to native units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- copy()
Return independent copy of this Reader.
New Reader will have its own file handle and can seek/iterate independently of the original.
Will also copy the current state of the Timestep held in the original Reader.
Changed in version 2.2.0: Arguments used to construct the reader are correctly captured and passed to the creation of the new class. Previously the only
n_atoms
was passed to class copies, leading to a class created with default parameters which may differ from the original class.
- property dimensions
unitcell dimensions (A, B, C, alpha, beta, gamma)
- property dt
timestep between frames
- property frame: int
Frame number of the current time step.
This is a simple short cut to
Timestep.frame
.
- get_aux_attribute(auxname, attrname)
Get the value of attrname from the auxiliary auxname
- Parameters:
See also
- get_aux_descriptions(auxnames=None)
Get descriptions to allow reloading the specified auxiliaries.
If no auxnames are provided, defaults to the full list of added auxiliaries.
Passing the resultant description to
add_auxiliary()
will allow recreation of the auxiliary. e.g., to duplicate all auxiliaries into a second trajectory:descriptions = trajectory_1.get_aux_descriptions() for aux in descriptions: trajectory_2.add_auxiliary(**aux)
- Returns:
List of dictionaries of the args/kwargs describing each auxiliary.
- Return type:
- iter_as_aux(auxname)
Iterate through timesteps for which there is at least one assigned step from the auxiliary auxname within the cutoff specified in auxname.
See also
- iter_auxiliary(auxname, start=None, stop=None, step=None, selected=None)
Iterate through the auxiliary auxname independently of the trajectory.
Will iterate over the specified steps of the auxiliary (defaults to all steps). Allows to access all values in an auxiliary, including those out of the time range of the trajectory, without having to also iterate through the trajectory.
After interation, the auxiliary will be repositioned at the current step.
- Parameters:
auxname (str) – Name of the auxiliary to iterate over.
(start (optional) – Options for iterating over a slice of the auxiliary.
stop (optional) – Options for iterating over a slice of the auxiliary.
step) (optional) – Options for iterating over a slice of the auxiliary.
selected (lst | ndarray, optional) – List of steps to iterate over.
- Yields:
AuxStep
object
See also
- property n_frames
number of frames in trajectory
- next_as_aux(auxname)
Move to the next timestep for which there is at least one step from the auxiliary auxname within the cutoff specified in auxname.
This allows progression through the trajectory without encountering
NaN
representative values (unless these are specifically part of the auxiliary data).If the auxiliary cutoff is not set, where auxiliary steps are less frequent (
auxiliary.dt > trajectory.dt
), this allows progression at the auxiliary pace (rounded to nearest timestep); while if the auxiliary steps are more frequent, this will work the same as callingnext()
.See the Auxiliary API.
See also
- static parse_n_atoms(filename, **kwargs)[source]
Read the coordinate file and deduce the number of atoms
- Returns:
n_atoms – the number of atoms in the coordinate file
- Return type:
- Raises:
NotImplementedError – when the number of atoms can’t be deduced
- rename_aux(auxname, new)
Change the name of the auxiliary auxname to new.
Provided there is not already an auxiliary named new, the auxiliary name will be changed in ts.aux namespace, the trajectory’s list of added auxiliaries, and in the auxiliary reader itself.
- Parameters:
- Raises:
ValueError – If the name new is already in use by an existing auxiliary.
- set_aux_attribute(auxname, attrname, new)
Set the value of attrname in the auxiliary auxname.
- Parameters:
See also
- property time
Time of the current frame in MDAnalysis time units (typically ps).
This is either read straight from the Timestep, or calculated as time =
Timestep.frame
*Timestep.dt
- timeseries(asel=None, atomgroup=None, start=None, stop=None, step=None, order='afc')[source]
Return a subset of coordinate data for an AtomGroup
- Parameters:
asel (
AtomGroup
) –The
AtomGroup
to read the coordinates from. Defaults to None, in which case the full set of coordinate data is returned.Deprecated since version 2.7.0: asel argument will be renamed to atomgroup in 3.0.0
atomgroup (AtomGroup (optional)) – Same as asel, will replace asel in 3.0.0
start (int (optional)) – Begin reading the trajectory at frame index start (where 0 is the index of the first frame in the trajectory); the default
None
starts at the beginning.stop (int (optional)) – End reading the trajectory at frame index stop-1, i.e, stop is excluded. The trajectory is read to the end with the default
None
.step (int (optional)) – Step size for reading; the default
None
is equivalent to 1 and means to read every frame.order (str (optional)) – the order/shape of the return data array, corresponding to (a)tom, (f)rame, (c)oordinates all six combinations of ‘a’, ‘f’, ‘c’ are allowed ie “fac” - return array where the shape is (frame, number of atoms, coordinates)
Changed in version 1.0.0: skip and format keywords have been removed.
Changed in version 2.4.0: ValueError now raised instead of NoDataError for empty input AtomGroup
- property totaltime: float
Total length of the trajectory
The time is calculated as
(n_frames - 1) * dt
, i.e., we assume that the first frame no time as elapsed. Thus, a trajectory with two frames will be considered to have a length of a single time step dt and a “trajectory” with a single frame will be reported as length 0.
- property transformations
Returns the list of transformations
- units = {'length': 'Angstrom', 'time': 'AKMA'}
dict with units of of time and length (and velocity, force, … for formats that support it)
- class MDAnalysis.coordinates.LAMMPS.DCDWriter(*args, **kwargs)[source]
Write a LAMMPS DCD trajectory.
The units can be set from the constructor with the keyword arguments timeunit and lengthunit. The defaults are “fs” and “Angstrom”. See
MDAnalysis.units
for other recognized values.- Parameters:
filename (str) – filename of trajectory
n_atoms (int) – number of atoms to be written
convert_units (bool (optional)) – convert from MDAnalysis units to format specific units
step (int (optional)) – number of steps between frames to be written
dt (float (optional)) – time between two frames. If
None
guess from first writtenTimeStep
remarks (str (optional)) – remarks to be stored in DCD. Shouldn’t be more then 240 characters
nsavc (int (optional)) – DCD files can also store the number of integrator time steps that correspond to the interval between two frames as nsavc (i.e., every how many MD steps is a frame saved to the DCD). By default, this number is just set to one and this should be sufficient for almost all cases but if required, nsavc can be changed.
istart (int (optional)) – starting frame number in integrator timesteps. CHARMM defaults to nsavc, i.e., start at frame number 1 = istart / nsavc. The value
None
will set istart to nsavc (the CHARMM default). The MDAnalysis default is 0 so that the frame number and time of the first frame is 0.**kwargs (dict) – General writer arguments
- convert_dimensions_to_unitcell(ts, inplace=True)
Read dimensions from timestep ts and return appropriate unitcell.
The default is to return
[A,B,C,alpha,beta,gamma]
; if this is not appropriate then this method has to be overriden.
- convert_forces_from_native(force, inplace=True)
Conversion of forces array force from native to base units
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_forces_to_native(force, inplace=True)
Conversion of force array force from base to native units.
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_pos_from_native(x, inplace=True)
Conversion of coordinate array x from native units to base units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_pos_to_native(x, inplace=True)
Conversion of coordinate array x from base units to native units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_from_native(t, inplace=True)
Convert time t from native units to base units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.) In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_to_native(t, inplace=True)
Convert time t from base units to native units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool, optional) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_velocities_from_native(v, inplace=True)
Conversion of velocities array v from native to base units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- convert_velocities_to_native(v, inplace=True)
Conversion of coordinate array v from base to native units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- has_valid_coordinates(criteria, x)
Returns
True
if all values are within limit values of their formats.Due to rounding, the test is asymmetric (and min is supposed to be negative):
min < x <= max
- Parameters:
criteria (dict) – dictionary containing the max and min values in native units
x (numpy.ndarray) –
(x, y, z)
coordinates of atoms selected to be written out
- Return type:
- units = {'length': 'Angstrom', 'time': 'AKMA'}
dict with units of of time and length (and velocity, force, … for formats that support it)
- write(obj)
Write current timestep, using the supplied obj.
Note
The size of the obj must be the same as the number of atoms provided when setting up the trajectory.
Changed in version 2.0.0: Deprecated support for Timestep argument to write has now been removed. Use AtomGroup or Universe as an input instead.
- class MDAnalysis.coordinates.LAMMPS.DATAReader(filename, **kwargs)[source]
Reads a single frame of coordinate information from a LAMMPS DATA file.
Added in version 0.9.0.
Changed in version 0.11.0: Frames now 0-based instead of 1-based
- OtherWriter(filename, **kwargs)
Returns a writer appropriate for filename.
Sets the default keywords start, step and dt (if available). n_atoms is always set from
Reader.n_atoms
.See also
Reader.Writer()
- Writer(filename, **kwargs)
A trajectory writer with the same properties as this trajectory.
- add_auxiliary(aux_spec: str | Dict[str, str] | None = None, auxdata: str | AuxReader | None = None, format: str | None = None, **kwargs) None
Add auxiliary data to be read alongside trajectory.
Auxiliary data may be any data timeseries from the trajectory additional to that read in by the trajectory reader. auxdata can be an
AuxReader
instance, or the data itself as e.g. a filename; in the latter case an appropriateAuxReader
is guessed from the data/file format. An appropriate format may also be directly provided as a key word argument.On adding, the AuxReader is initially matched to the current timestep of the trajectory, and will be updated when the trajectory timestep changes (through a call to
next()
or jumping timesteps withtrajectory[i]
).The representative value(s) of the auxiliary data for each timestep (as calculated by the
AuxReader
) are stored in the current timestep in thets.aux
namespace under aux_spec; e.g. to add additional pull force data stored in pull-force.xvg:u = MDAnalysis.Universe(PDB, XTC) u.trajectory.add_auxiliary('pull', 'pull-force.xvg')
The representative value for the current timestep may then be accessed as
u.trajectory.ts.aux.pull
oru.trajectory.ts.aux['pull']
.The following applies to energy readers like the
EDRReader
.All data that is present in the (energy) file can be added by omitting aux_spec like so:
u.trajectory.add_auxiliary(auxdata="ener.edr")
aux_spec is expected to be a dictionary that maps the desired attribute name in the
ts.aux
namespace to the precise data to be added as identified by adata_selector
:term_dict = {"temp": "Temperature", "epot": "Potential"} u.trajectory.add_auxiliary(term_dict, "ener.edr")
Adding this data can be useful, for example, to filter trajectory frames based on non-coordinate data like the potential energy of each time step. Trajectory slicing allows working on a subset of frames:
selected_frames = np.array([ts.frame for ts in u.trajectory if ts.aux.epot < some_threshold]) subset = u.trajectory[selected_frames]
See also
Note
Auxiliary data is assumed to be time-ordered, with no duplicates. See the Auxiliary API.
- add_transformations(*transformations)
Add all transformations to be applied to the trajectory.
This function take as list of transformations as an argument. These transformations are functions that will be called by the Reader and given a
Timestep
object as argument, which will be transformed and returned to the Reader. The transformations can be part of thetransformations
module, or created by the user, and are stored as a list transformations. This list can only be modified once, and further calls of this function will raise an exception.u = MDAnalysis.Universe(topology, coordinates) workflow = [some_transform, another_transform, this_transform] u.trajectory.add_transformations(*workflow)
- Parameters:
transform_list (list) – list of all the transformations that will be applied to the coordinates
See also
- property aux_list
Lists the names of added auxiliary data.
- check_slice_indices(start, stop, step)
Check frame indices are valid and clip to fit trajectory.
The usage follows standard Python conventions for
range()
but see the warning below.- Parameters:
start (int or None) – Starting frame index (inclusive).
None
corresponds to the default of 0, i.e., the initial frame.stop (int or None) – Last frame index (exclusive).
None
corresponds to the default of n_frames, i.e., it includes the last frame of the trajectory.step (int or None) – step size of the slice,
None
corresponds to the default of 1, i.e, include every frame in the range start, stop.
- Returns:
start, stop, step – Integers representing the slice
- Return type:
Warning
The returned values start, stop and step give the expected result when passed in
range()
but gives unexpected behavior when passed in aslice
whenstop=None
andstep=-1
This can be a problem for downstream processing of the output from this method. For example, slicing of trajectories is implemented by passing the values returned by
check_slice_indices()
torange()
range(start, stop, step)
and using them as the indices to randomly seek to. On the other hand, in
MDAnalysis.analysis.base.AnalysisBase
the values returned bycheck_slice_indices()
are used to splice the trajectory by creating aslice
instanceslice(start, stop, step)
This creates a discrepancy because these two lines are not equivalent:
range(10, -1, -1) # [10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0] range(10)[slice(10, -1, -1)] # []
- close()
Close the trajectory file.
- convert_forces_from_native(force, inplace=True)
Conversion of forces array force from native to base units
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_forces_to_native(force, inplace=True)
Conversion of force array force from base to native units.
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_pos_from_native(x, inplace=True)
Conversion of coordinate array x from native units to base units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_pos_to_native(x, inplace=True)
Conversion of coordinate array x from base units to native units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_from_native(t, inplace=True)
Convert time t from native units to base units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.) In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_to_native(t, inplace=True)
Convert time t from base units to native units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool, optional) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_velocities_from_native(v, inplace=True)
Conversion of velocities array v from native to base units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- convert_velocities_to_native(v, inplace=True)
Conversion of coordinate array v from base to native units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- copy()
Return independent copy of this Reader.
New Reader will have its own file handle and can seek/iterate independently of the original.
Will also copy the current state of the Timestep held in the original Reader.
Changed in version 2.2.0: Arguments used to construct the reader are correctly captured and passed to the creation of the new class. Previously the only
n_atoms
was passed to class copies, leading to a class created with default parameters which may differ from the original class.
- property frame: int
Frame number of the current time step.
This is a simple short cut to
Timestep.frame
.
- get_aux_attribute(auxname, attrname)
Get the value of attrname from the auxiliary auxname
- Parameters:
See also
- get_aux_descriptions(auxnames=None)
Get descriptions to allow reloading the specified auxiliaries.
If no auxnames are provided, defaults to the full list of added auxiliaries.
Passing the resultant description to
add_auxiliary()
will allow recreation of the auxiliary. e.g., to duplicate all auxiliaries into a second trajectory:descriptions = trajectory_1.get_aux_descriptions() for aux in descriptions: trajectory_2.add_auxiliary(**aux)
- Returns:
List of dictionaries of the args/kwargs describing each auxiliary.
- Return type:
- iter_as_aux(auxname)
Iterate through timesteps for which there is at least one assigned step from the auxiliary auxname within the cutoff specified in auxname.
See also
- iter_auxiliary(auxname, start=None, stop=None, step=None, selected=None)
Iterate through the auxiliary auxname independently of the trajectory.
Will iterate over the specified steps of the auxiliary (defaults to all steps). Allows to access all values in an auxiliary, including those out of the time range of the trajectory, without having to also iterate through the trajectory.
After interation, the auxiliary will be repositioned at the current step.
- Parameters:
auxname (str) – Name of the auxiliary to iterate over.
(start (optional) – Options for iterating over a slice of the auxiliary.
stop (optional) – Options for iterating over a slice of the auxiliary.
step) (optional) – Options for iterating over a slice of the auxiliary.
selected (lst | ndarray, optional) – List of steps to iterate over.
- Yields:
AuxStep
object
See also
- next()
Forward one step to next frame.
- next_as_aux(auxname)
Move to the next timestep for which there is at least one step from the auxiliary auxname within the cutoff specified in auxname.
This allows progression through the trajectory without encountering
NaN
representative values (unless these are specifically part of the auxiliary data).If the auxiliary cutoff is not set, where auxiliary steps are less frequent (
auxiliary.dt > trajectory.dt
), this allows progression at the auxiliary pace (rounded to nearest timestep); while if the auxiliary steps are more frequent, this will work the same as callingnext()
.See the Auxiliary API.
See also
- classmethod parse_n_atoms(filename, **kwargs)
Read the coordinate file and deduce the number of atoms
- Returns:
n_atoms – the number of atoms in the coordinate file
- Return type:
- Raises:
NotImplementedError – when the number of atoms can’t be deduced
- rename_aux(auxname, new)
Change the name of the auxiliary auxname to new.
Provided there is not already an auxiliary named new, the auxiliary name will be changed in ts.aux namespace, the trajectory’s list of added auxiliaries, and in the auxiliary reader itself.
- Parameters:
- Raises:
ValueError – If the name new is already in use by an existing auxiliary.
- rewind()
Position at beginning of trajectory
- set_aux_attribute(auxname, attrname, new)
Set the value of attrname in the auxiliary auxname.
- Parameters:
See also
- property time
Time of the current frame in MDAnalysis time units (typically ps).
This is either read straight from the Timestep, or calculated as time =
Timestep.frame
*Timestep.dt
- timeseries(asel: AtomGroup | None = None, atomgroup: Atomgroup | None = None, start: int | None = None, stop: int | None = None, step: int | None = None, order: str | None = 'fac') ndarray
Return a subset of coordinate data for an AtomGroup
- Parameters:
asel (AtomGroup (optional)) –
The
AtomGroup
to read the coordinates from. Defaults toNone
, in which case the full set of coordinate data is returned.Deprecated since version 2.7.0: asel argument will be renamed to atomgroup in 3.0.0
atomgroup (AtomGroup (optional)) – Same as asel, will replace asel in 3.0.0
start (int (optional)) – Begin reading the trajectory at frame index start (where 0 is the index of the first frame in the trajectory); the default
None
starts at the beginning.stop (int (optional)) – End reading the trajectory at frame index stop-1, i.e, stop is excluded. The trajectory is read to the end with the default
None
.step (int (optional)) – Step size for reading; the default
None
is equivalent to 1 and means to read every frame.order (str (optional)) – the order/shape of the return data array, corresponding to (a)tom, (f)rame, (c)oordinates all six combinations of ‘a’, ‘f’, ‘c’ are allowed ie “fac” - return array where the shape is (frame, number of atoms, coordinates)
See also
Added in version 2.4.0.
- property totaltime: float
Total length of the trajectory
The time is calculated as
(n_frames - 1) * dt
, i.e., we assume that the first frame no time as elapsed. Thus, a trajectory with two frames will be considered to have a length of a single time step dt and a “trajectory” with a single frame will be reported as length 0.
- property transformations
Returns the list of transformations
- units = {'length': 'Angstrom', 'time': None, 'velocity': 'Angstrom/fs'}
dict with units of of time and length (and velocity, force, … for formats that support it)
- class MDAnalysis.coordinates.LAMMPS.DATAWriter(filename, convert_units=True, **kwargs)[source]
Write out the current time step as a LAMMPS DATA file.
This writer supports the sections Atoms, Masses, Velocities, Bonds, Angles, Dihedrals, and Impropers. This writer will write the header and these sections (if applicable). Atoms section is written in the “full” sub-style if charges are available or “molecular” sub-style if they are not. Molecule id is set to 0 for all atoms.
Note
This writer assumes “conventional” or “real” LAMMPS units where length is measured in Angstroms and velocity is measured in Angstroms per femtosecond. To write in different units, specify lengthunit
If atom types are not already positive integers, the user must set them to be positive integers, because the writer will not automatically assign new types.
To preserve numerical atom types when writing a selection, the Masses section will have entries for each atom type up to the maximum atom type. If the universe does not contain atoms of some type in {1, … max(atom_types)}, then the mass for that type will be set to 1.
In order to write bonds, each selected bond type must be explicitly set to an integer >= 1.
Set up a DATAWriter
- Parameters:
- close()
Close the trajectory file.
- convert_dimensions_to_unitcell(ts, inplace=True)
Read dimensions from timestep ts and return appropriate unitcell.
The default is to return
[A,B,C,alpha,beta,gamma]
; if this is not appropriate then this method has to be overriden.
- convert_forces_from_native(force, inplace=True)
Conversion of forces array force from native to base units
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_forces_to_native(force, inplace=True)
Conversion of force array force from base to native units.
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_pos_from_native(x, inplace=True)
Conversion of coordinate array x from native units to base units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_pos_to_native(x, inplace=True)
Conversion of coordinate array x from base units to native units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_from_native(t, inplace=True)
Convert time t from native units to base units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.) In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_to_native(t, inplace=True)
Convert time t from base units to native units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool, optional) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_velocities_from_native(v, inplace=True)
Conversion of velocities array v from native to base units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- convert_velocities_to_native(v, inplace=True)
Conversion of coordinate array v from base to native units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- has_valid_coordinates(criteria, x)
Returns
True
if all values are within limit values of their formats.Due to rounding, the test is asymmetric (and min is supposed to be negative):
min < x <= max
- Parameters:
criteria (dict) – dictionary containing the max and min values in native units
x (numpy.ndarray) –
(x, y, z)
coordinates of atoms selected to be written out
- Return type:
- units = {'length': None, 'time': None, 'velocity': None}
dict with units of of time and length (and velocity, force, … for formats that support it)
- write(selection, frame=None)[source]
Write selection at current trajectory frame to file.
The sections for Atoms, Masses, Velocities, Bonds, Angles, Dihedrals, and Impropers (if these are defined) are written. The Atoms section is written in the “full” sub-style if charges are available or “molecular” sub-style if they are not. Molecule id in atoms section is set to to 0.
No other sections are written to the DATA file. As of this writing, other sections are not parsed into the topology by the
DATAReader
.Note
If the selection includes a partial fragment, then only the bonds, angles, etc. whose atoms are contained within the selection will be included.
- class MDAnalysis.coordinates.LAMMPS.DumpReader(filename, lammps_coordinate_convention='auto', unwrap_images=False, additional_columns=None, **kwargs)[source]
Reads the default LAMMPS dump format
Supports coordinates in the LAMMPS “unscaled” (x,y,z), “scaled” (xs,ys,zs), “unwrapped” (xu,yu,zu) and “scaled_unwrapped” (xsu,ysu,zsu) coordinate conventions (see https://docs.lammps.org/dump.html for more details). If lammps_coordinate_convention=’auto’ (default), one will be guessed. Guessing checks whether the coordinates fit each convention in the order “unscaled”, “scaled”, “unwrapped”, “scaled_unwrapped” and whichever set of coordinates is detected first will be used. If coordinates are given in the scaled coordinate convention (xs,ys,zs) or scaled unwrapped coordinate convention (xsu,ysu,zsu) they will automatically be converted from their scaled/fractional representation to their real values.
Supports both orthogonal and triclinic simulation box dimensions (for more details see https://docs.lammps.org/Howto_triclinic.html). In either case, MDAnalysis will always use
(*A*, *B*, *C*, *alpha*, *beta*, *gamma*)
to represent the unit cell. Lengths A, B, C are in the MDAnalysis length unit (Å), and angles are in degrees.By using the keyword additional_columns, you can specify arbitrary data to be read. The keyword expects a list of the names of the columns or True to read all additional columns. The results are saved to
Timestep.data
. For example, if your LAMMPS dump looks like thisITEM: ATOMS id x y z q l 1 2.84 8.17 -25 0.00258855 1.1 2 7.1 8.17 -25 6.91952e-05 1.2
Then you may parse the additional columns q and l via:
u = mda.Universe('structure.data', 'traj.lammpsdump', additional_columns=['q', 'l'])
The additional data is then available for each time step via:
for ts in u.trajectory: charges = ts.data['q'] # Access additional data, sorted by the id ls = ts.data['l'] ...
- Parameters:
filename (str) – Filename of LAMMPS dump file
lammps_coordinate_convention (str (optional) default="auto") –
Convention used in coordinates, can be one of the following according to the LAMMPS documentation:
”auto” - Detect coordinate type from file column header. If auto detection is used, the guessing checks whether the coordinates fit each convention in the order “unscaled”, “scaled”, “unwrapped”, “scaled_unwrapped” and whichever set of coordinates is detected first will be used.
- ”scaled” - Coordinates wrapped in box and scaled by box length (see
note below), i.e., xs, ys, zs
”scaled_unwrapped” - Coordinates unwrapped and scaled by box length, (see note below) i.e., xsu, ysu, zsu
”unscaled” - Coordinates wrapped in box, i.e., x, y, z
”unwrapped” - Coordinates unwrapped, i.e., xu, yu, zu
If coordinates are given in the scaled coordinate convention (xs,ys,zs) or scaled unwrapped coordinate convention (xsu,ysu,zsu) they will automatically be converted from their scaled/fractional representation to their real values.
unwrap_images (bool (optional) default=False) – If True and the dump file contains image flags, the coordinates will be unwrapped. See read_data in the lammps documentation for more information.
**kwargs – Other keyword arguments used in
ReaderBase
Changed in version 2.7.0: Reading of arbitrary, additional columns is now supported. (Issue #3608)
Changed in version 2.4.0: Now imports velocities and forces, translates the box to the origin, and optionally unwraps trajectories with image flags upon loading.
Changed in version 2.2.0: Triclinic simulation boxes are supported. (Issue #3383)
Changed in version 2.0.0: Now parses coordinates in multiple lammps conventions (x,xs,xu,xsu)
Added in version 0.19.0.
- OtherWriter(filename, **kwargs)
Returns a writer appropriate for filename.
Sets the default keywords start, step and dt (if available). n_atoms is always set from
Reader.n_atoms
.See also
Reader.Writer()
- Writer(filename, **kwargs)
A trajectory writer with the same properties as this trajectory.
- add_auxiliary(aux_spec: str | Dict[str, str] | None = None, auxdata: str | AuxReader | None = None, format: str | None = None, **kwargs) None
Add auxiliary data to be read alongside trajectory.
Auxiliary data may be any data timeseries from the trajectory additional to that read in by the trajectory reader. auxdata can be an
AuxReader
instance, or the data itself as e.g. a filename; in the latter case an appropriateAuxReader
is guessed from the data/file format. An appropriate format may also be directly provided as a key word argument.On adding, the AuxReader is initially matched to the current timestep of the trajectory, and will be updated when the trajectory timestep changes (through a call to
next()
or jumping timesteps withtrajectory[i]
).The representative value(s) of the auxiliary data for each timestep (as calculated by the
AuxReader
) are stored in the current timestep in thets.aux
namespace under aux_spec; e.g. to add additional pull force data stored in pull-force.xvg:u = MDAnalysis.Universe(PDB, XTC) u.trajectory.add_auxiliary('pull', 'pull-force.xvg')
The representative value for the current timestep may then be accessed as
u.trajectory.ts.aux.pull
oru.trajectory.ts.aux['pull']
.The following applies to energy readers like the
EDRReader
.All data that is present in the (energy) file can be added by omitting aux_spec like so:
u.trajectory.add_auxiliary(auxdata="ener.edr")
aux_spec is expected to be a dictionary that maps the desired attribute name in the
ts.aux
namespace to the precise data to be added as identified by adata_selector
:term_dict = {"temp": "Temperature", "epot": "Potential"} u.trajectory.add_auxiliary(term_dict, "ener.edr")
Adding this data can be useful, for example, to filter trajectory frames based on non-coordinate data like the potential energy of each time step. Trajectory slicing allows working on a subset of frames:
selected_frames = np.array([ts.frame for ts in u.trajectory if ts.aux.epot < some_threshold]) subset = u.trajectory[selected_frames]
See also
Note
Auxiliary data is assumed to be time-ordered, with no duplicates. See the Auxiliary API.
- add_transformations(*transformations)
Add all transformations to be applied to the trajectory.
This function take as list of transformations as an argument. These transformations are functions that will be called by the Reader and given a
Timestep
object as argument, which will be transformed and returned to the Reader. The transformations can be part of thetransformations
module, or created by the user, and are stored as a list transformations. This list can only be modified once, and further calls of this function will raise an exception.u = MDAnalysis.Universe(topology, coordinates) workflow = [some_transform, another_transform, this_transform] u.trajectory.add_transformations(*workflow)
The transformations are applied in the order given in the list transformations, i.e., the first transformation is the first or innermost one to be applied to the
Timestep
. The example above would be equivalent tofor ts in u.trajectory: ts = this_transform(another_transform(some_transform(ts)))
- Parameters:
transform_list (list) – list of all the transformations that will be applied to the coordinates in the order given in the list
See also
- property aux_list
Lists the names of added auxiliary data.
- check_slice_indices(start, stop, step)
Check frame indices are valid and clip to fit trajectory.
The usage follows standard Python conventions for
range()
but see the warning below.- Parameters:
start (int or None) – Starting frame index (inclusive).
None
corresponds to the default of 0, i.e., the initial frame.stop (int or None) – Last frame index (exclusive).
None
corresponds to the default of n_frames, i.e., it includes the last frame of the trajectory.step (int or None) – step size of the slice,
None
corresponds to the default of 1, i.e, include every frame in the range start, stop.
- Returns:
start, stop, step – Integers representing the slice
- Return type:
Warning
The returned values start, stop and step give the expected result when passed in
range()
but gives unexpected behavior when passed in aslice
whenstop=None
andstep=-1
This can be a problem for downstream processing of the output from this method. For example, slicing of trajectories is implemented by passing the values returned by
check_slice_indices()
torange()
range(start, stop, step)
and using them as the indices to randomly seek to. On the other hand, in
MDAnalysis.analysis.base.AnalysisBase
the values returned bycheck_slice_indices()
are used to splice the trajectory by creating aslice
instanceslice(start, stop, step)
This creates a discrepancy because these two lines are not equivalent:
range(10, -1, -1) # [10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0] range(10)[slice(10, -1, -1)] # []
- convert_forces_from_native(force, inplace=True)
Conversion of forces array force from native to base units
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_forces_to_native(force, inplace=True)
Conversion of force array force from base to native units.
- Parameters:
force (array_like) – Forces to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.7.
- convert_pos_from_native(x, inplace=True)
Conversion of coordinate array x from native units to base units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_pos_to_native(x, inplace=True)
Conversion of coordinate array x from base units to native units.
- Parameters:
x (array_like) – Positions to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_from_native(t, inplace=True)
Convert time t from native units to base units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.) In-place operations improve performance because allocating new arrays is avoided.
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_time_to_native(t, inplace=True)
Convert time t from base units to native units.
- Parameters:
t (array_like) – Time values to transform
inplace (bool, optional) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to
False
so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
- convert_velocities_from_native(v, inplace=True)
Conversion of velocities array v from native to base units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- convert_velocities_to_native(v, inplace=True)
Conversion of coordinate array v from base to native units
- Parameters:
v (array_like) – Velocities to transform
inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
Note
By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.
Added in version 0.7.5.
- copy()
Return independent copy of this Reader.
New Reader will have its own file handle and can seek/iterate independently of the original.
Will also copy the current state of the Timestep held in the original Reader.
Changed in version 2.2.0: Arguments used to construct the reader are correctly captured and passed to the creation of the new class. Previously the only
n_atoms
was passed to class copies, leading to a class created with default parameters which may differ from the original class.
- property frame: int
Frame number of the current time step.
This is a simple short cut to
Timestep.frame
.
- get_aux_attribute(auxname, attrname)
Get the value of attrname from the auxiliary auxname
- Parameters:
See also
- get_aux_descriptions(auxnames=None)
Get descriptions to allow reloading the specified auxiliaries.
If no auxnames are provided, defaults to the full list of added auxiliaries.
Passing the resultant description to
add_auxiliary()
will allow recreation of the auxiliary. e.g., to duplicate all auxiliaries into a second trajectory:descriptions = trajectory_1.get_aux_descriptions() for aux in descriptions: trajectory_2.add_auxiliary(**aux)
- Returns:
List of dictionaries of the args/kwargs describing each auxiliary.
- Return type:
- iter_as_aux(auxname)
Iterate through timesteps for which there is at least one assigned step from the auxiliary auxname within the cutoff specified in auxname.
See also
- iter_auxiliary(auxname, start=None, stop=None, step=None, selected=None)
Iterate through the auxiliary auxname independently of the trajectory.
Will iterate over the specified steps of the auxiliary (defaults to all steps). Allows to access all values in an auxiliary, including those out of the time range of the trajectory, without having to also iterate through the trajectory.
After interation, the auxiliary will be repositioned at the current step.
- Parameters:
auxname (str) – Name of the auxiliary to iterate over.
(start (optional) – Options for iterating over a slice of the auxiliary.
stop (optional) – Options for iterating over a slice of the auxiliary.
step) (optional) – Options for iterating over a slice of the auxiliary.
selected (lst | ndarray, optional) – List of steps to iterate over.
- Yields:
AuxStep
object
See also
- next_as_aux(auxname)
Move to the next timestep for which there is at least one step from the auxiliary auxname within the cutoff specified in auxname.
This allows progression through the trajectory without encountering
NaN
representative values (unless these are specifically part of the auxiliary data).If the auxiliary cutoff is not set, where auxiliary steps are less frequent (
auxiliary.dt > trajectory.dt
), this allows progression at the auxiliary pace (rounded to nearest timestep); while if the auxiliary steps are more frequent, this will work the same as callingnext()
.See the Auxiliary API.
See also
- classmethod parse_n_atoms(filename, **kwargs)
Read the coordinate file and deduce the number of atoms
- Returns:
n_atoms – the number of atoms in the coordinate file
- Return type:
- Raises:
NotImplementedError – when the number of atoms can’t be deduced
- rename_aux(auxname, new)
Change the name of the auxiliary auxname to new.
Provided there is not already an auxiliary named new, the auxiliary name will be changed in ts.aux namespace, the trajectory’s list of added auxiliaries, and in the auxiliary reader itself.
- Parameters:
- Raises:
ValueError – If the name new is already in use by an existing auxiliary.
- set_aux_attribute(auxname, attrname, new)
Set the value of attrname in the auxiliary auxname.
- Parameters:
See also
- property time
Time of the current frame in MDAnalysis time units (typically ps).
This is either read straight from the Timestep, or calculated as time =
Timestep.frame
*Timestep.dt
- timeseries(asel: AtomGroup | None = None, atomgroup: Atomgroup | None = None, start: int | None = None, stop: int | None = None, step: int | None = None, order: str | None = 'fac') ndarray
Return a subset of coordinate data for an AtomGroup
- Parameters:
asel (AtomGroup (optional)) –
The
AtomGroup
to read the coordinates from. Defaults toNone
, in which case the full set of coordinate data is returned.Deprecated since version 2.7.0: asel argument will be renamed to atomgroup in 3.0.0
atomgroup (AtomGroup (optional)) – Same as asel, will replace asel in 3.0.0
start (int (optional)) – Begin reading the trajectory at frame index start (where 0 is the index of the first frame in the trajectory); the default
None
starts at the beginning.stop (int (optional)) – End reading the trajectory at frame index stop-1, i.e, stop is excluded. The trajectory is read to the end with the default
None
.step (int (optional)) – Step size for reading; the default
None
is equivalent to 1 and means to read every frame.order (str (optional)) – the order/shape of the return data array, corresponding to (a)tom, (f)rame, (c)oordinates all six combinations of ‘a’, ‘f’, ‘c’ are allowed ie “fac” - return array where the shape is (frame, number of atoms, coordinates)
See also
Added in version 2.4.0.
- property totaltime: float
Total length of the trajectory
The time is calculated as
(n_frames - 1) * dt
, i.e., we assume that the first frame no time as elapsed. Thus, a trajectory with two frames will be considered to have a length of a single time step dt and a “trajectory” with a single frame will be reported as length 0.
- property transformations
Returns the list of transformations
- units = {'length': None, 'time': None, 'velocity': None}
dict with units of of time and length (and velocity, force, … for formats that support it)