7.33. Base classes — MDAnalysis.coordinates.base
Derive, FrameIterator, Reader and Writer classes from the classes in this module. The derived classes must follow the Trajectory API.
7.33.1. FrameIterators
FrameIterators are “sliced trajectories” (a trajectory is a Reader) that can be iterated over. They are typically created by slicing a trajectory or by fancy-indexing of a trajectory with an array of frame numbers or a boolean mask of all frames.
Iterator classes used by the by the ProtoReader
:
- class MDAnalysis.coordinates.base.FrameIteratorBase(trajectory)[source]
Base iterable over the frames of a trajectory.
A frame iterable has a length that can be accessed with the
len()
function, and can be indexed similarly to a full trajectory. When indexed, indices are resolved relative to the iterable and not relative to the trajectory.Added in version 0.19.0.
- class MDAnalysis.coordinates.base.FrameIteratorSliced(trajectory, frames)[source]
Iterable over the frames of a trajectory on the basis of a slice.
- Parameters:
trajectory (ProtoReader) – The trajectory over which to iterate.
frames (slice) – A slice to select the frames of interest.
See also
- class MDAnalysis.coordinates.base.FrameIteratorAll(trajectory)[source]
Iterable over all the frames of a trajectory.
- Parameters:
trajectory (ProtoReader) – The trajectory over which to iterate.
See also
- class MDAnalysis.coordinates.base.FrameIteratorIndices(trajectory, frames)[source]
Iterable over the frames of a trajectory listed in a sequence of indices.
- Parameters:
trajectory (ProtoReader) – The trajectory over which to iterate.
frames (sequence) – A sequence of indices.
See also
- class MDAnalysis.coordinates.base.StreamFrameIteratorSliced(trajectory, step)[source]
Iterator for sliced frames in a streamed trajectory.
Created when slicing a streaming trajectory with a step parameter (e.g.,
trajectory[::n]
). Reads every nth frame from the continuous stream, discarding intermediate frames for performance.This differs from iterating over all frames (
trajectory[:]
) which usesFrameIteratorAll
and processes every frame sequentially without skipping.Streaming constraints apply to the sliced iterator:
Frames cannot be accessed randomly (no indexing support)
The total number of frames is unknown until streaming ends
Rewinding or restarting iteration is not possible
Only forward iteration with a fixed step size is supported
- Parameters:
trajectory (StreamReaderBase) – The streaming trajectory reader to iterate over. Must be a stream-based reader that supports continuous data reading.
step (int) – Step size for iteration. Must be a positive integer. A step of 1 reads every frame, step of 2 reads every other frame, etc.
See also
Added in version 2.10.0.
7.33.2. Readers
Readers know how to take trajectory data in a given format and present it in a common API to the user in MDAnalysis. There are two types of readers:
Readers for multi frame trajectories, i.e., file formats that typically contain many frames. These readers are typically derived from
ReaderBase
.Readers for single frame formats: These file formats only contain a single coordinate set. These readers are derived from
SingleFrameReaderBase
.
The underlying low-level readers handle closing of files in different
ways. Typically, the MDAnalysis readers try to ensure that files are always
closed when a reader instance is garbage collected, which relies on
implementing a __del__()
method. However, in some cases, this
is not necessary (for instance, for the single frame formats) and then such a
method can lead to undesirable side effects (such as memory leaks). In this
case, ProtoReader
should be used.
- class MDAnalysis.coordinates.base.ReaderBase(filename, convert_units=True, **kwargs)[source]
Base class for trajectory readers that extends
ProtoReader
with a__del__()
method.New Readers should subclass
ReaderBase
and properly implement aclose()
method, to ensure proper release of resources (mainly file handles). Readers that are inherently safe in this regard should subclassProtoReader
instead.See the Trajectory API definition in for the required attributes and methods.
See also
Changed in version 0.11.0: Most of the base Reader class definitions were offloaded to
ProtoReader
so as to allow the subclassing of ReaderBases without a__del__()
method. Created init method to create common functionality, all ReaderBase subclasses must nowsuper()
through this class. Added attribute_ts_kwargs
, which is created in init. Provides kwargs to be passed toTimestep
Changed in version 1.0: Removed deprecated flags functionality, use convert_units kwarg instead
- 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)] # []
- 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()[source]
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)
- class MDAnalysis.coordinates.base.SingleFrameReaderBase(filename, convert_units=True, n_atoms=None, **kwargs)[source]
Base class for Readers that only have one frame.
To use this base class, define the method
_read_first_frame()
to read from file self.filename. This should populate the attribute self.ts with aTimestep
object.Added in version 0.10.0.
Changed in version 0.11.0: Added attribute “_ts_kwargs” for subclasses Keywords “dt” and “time_offset” read into _ts_kwargs
Changed in version 2.2.0: Calling __iter__ now rewinds the reader before yielding a
Timestep
object (fixing behavior that was not well defined previously).Changed in version 2.10.0: Fixed a typo in the attribute assignment (self.atom → self.atoms), which may affect subclasses relying on this value.
- 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)[source]
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)] # []
- 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()[source]
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)
- class MDAnalysis.coordinates.base.ProtoReader[source]
Base class for Readers, without a
__del__()
method.Extends
IOBase
with most attributes and methods of a generic Reader, with the exception of a__del__()
method. It should be used as base for Readers that do not need__del__()
, especially since having even an empty__del__()
might lead to memory leaks.See the Trajectory API definition in
MDAnalysis.coordinates.__init__
for the required attributes and methods.See also
Changed in version 0.11.0: Frames now 0-based instead of 1-based
Changed in version 2.0.0: Now supports (un)pickle. Upon unpickling, the current timestep is retained by reconstrunction.
Changed in version 2.8.0: the modification of coordinates was preserved after serialization.
- OtherWriter(filename, **kwargs)[source]
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)[source]
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 [source]
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)[source]
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)[source]
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)] # []
- property frame: int
Frame number of the current time step.
This is a simple short cut to
Timestep.frame
.
- get_aux_attribute(auxname, attrname)[source]
Get the value of attrname from the auxiliary auxname
- Parameters:
See also
- get_aux_descriptions(auxnames=None)[source]
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)[source]
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)[source]
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)[source]
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)[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
- remove_auxiliary(auxname)[source]
Clear data and close the
AuxReader
for the auxiliary auxname.See also
- rename_aux(auxname, new)[source]
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)[source]
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 [source]
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
- class MDAnalysis.coordinates.base.StreamReaderBase(filename, convert_units=True, **kwargs)[source]
Base class for readers that read a continuous stream of data.
This class is designed for readers that process continuous data streams, such as live feeds from simulations. Unlike traditional trajectory readers that can randomly access frames, streaming readers have fundamental constraints:
No random access: Cannot seek to arbitrary frames (no
traj[5]
)Forward-only: Can only iterate sequentially through frames
No length: Total number of frames is unknown until stream ends
No rewinding: Cannot restart or rewind the stream
No copying: Cannot create independent copies of the reader
No reopening: Cannot restart iteration once stream is consumed
No timeseries: Cannot use
timeseries()
or bulk data extractionNo writers: Cannot create
Writer()
orOtherWriter()
instancesNo pickling: Cannot serialize reader instances (limits multiprocessing)
No StreamWriterBase: No complementary Writer class available for streaming data
The reader raises
RuntimeError
for operations that require random access or rewinding, includingrewind()
,copy()
,timeseries()
,Writer()
,OtherWriter()
, andlen()
. Only slice notation is supported for iteration.- Parameters:
See also
StreamFrameIteratorSliced
Iterator for stepped streaming access
ReaderBase
Base class for standard trajectory readers
Added in version 2.10.0.
- OtherWriter(filename, **kwargs)[source]
Writer creation is not supported for streaming trajectories.
OtherWriter initialization requires frame-based parameters and trajectory indexing information. Streaming readers process data sequentially without meaningful frame indexing, making writer setup impossible.
- Parameters:
filename (str) – Output filename (ignored, as method is not supported)
**kwargs – Additional keyword arguments (ignored, as method is not supported)
- Raises:
RuntimeError – Always raised, as writer creation is not supported for streaming trajectories
- Writer(filename, **kwargs)[source]
Writer creation is not supported for streaming trajectories.
Writer creation requires trajectory metadata that streaming readers cannot provide due to their sequential processing nature.
- Parameters:
filename (str) – Output filename (ignored, as method is not supported)
**kwargs – Additional keyword arguments (ignored, as method is not supported)
- Raises:
RuntimeError – Always raised, as writer creation is not supported for streaming trajectories
- check_slice_indices(start, stop, step)[source]
Check and validate slice indices for streaming trajectories.
Streaming trajectories have fundamental constraints that differ from traditional trajectory files:
No start/stop indices: Since streams process data continuously without knowing the total length,
start
andstop
must beNone
Step-only slicing: Only the
step
parameter is meaningful, controlling how many frames to skip during iterationForward-only:
step
must be positive (> 0) as streams cannot be processed backward in time
- Parameters:
- Returns:
(start, stop, step) with validated values
- Return type:
- Raises:
ValueError – If
start
orstop
are notNone
, or ifstep
is not a positive integer.
Examples
Valid streaming slices:
traj[:] # All frames (step=None, equivalent to step=1) traj[::2] # Every 2nd frame traj[::10] # Every 10th frame
Invalid streaming slices:
traj[5:] # Cannot specify start index traj[:100] # Cannot specify stop index traj[5:100:2] # Cannot specify start or stop indices traj[::-1] # Cannot go backwards (negative step)
See also
__getitem__
,StreamFrameIteratorSliced
Added in version 2.10.0.
- copy()[source]
Reader copying is not supported for streaming trajectories.
Streaming readers maintain internal state and connection resources that cannot be duplicated. Each stream connection is unique and cannot be copied.
- Raises:
RuntimeError – Always raised, as copying is not supported for streaming trajectories
- property n_frames
Changes as stream is processed unlike other readers
- next()[source]
Advance to the next timestep in the streaming trajectory.
Streaming readers process frames sequentially and cannot rewind once iteration completes. Use
for ts in trajectory
for iteration.- Returns:
The next timestep in the stream
- Return type:
- Raises:
StopIteration – When the stream ends or no more frames are available
- rewind()[source]
Rewinding is not supported for streaming trajectories.
Streaming readers process data continuously from streams and cannot restart or go backward in the stream once consumed.
- Raises:
RuntimeError – Always raised, as rewinding is not supported for streaming trajectories
- timeseries(**kwargs)[source]
Timeseries extraction is not supported for streaming trajectories.
Streaming readers cannot randomly access frames or store bulk coordinate data in memory, which
timeseries()
requires. Use sequential frame iteration instead.- Parameters:
**kwargs – Any keyword arguments (ignored, as method is not supported)
- Raises:
RuntimeError – Always raised, as timeseries extraction is not supported for streaming trajectories
7.33.3. Writers
Writers know how to write information in a Timestep
to a trajectory
file.
- class MDAnalysis.coordinates.base.WriterBase[source]
Base class for trajectory writers.
See Trajectory API definition in for the required attributes and methods.
Changed in version 2.0.0: Deprecated
write_next_timestep()
has now been removed, please usewrite()
instead.- close()
Close the trajectory file.
- convert_dimensions_to_unitcell(ts, inplace=True)[source]
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)[source]
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(obj)[source]
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.
7.33.4. Converters
Converters output information to other libraries.
Deprecated since version 2.7.0: All converter code has been moved to MDAnalysis.converters
and will
be removed from the MDAnalysis.coordinates.base
module in 3.0.0.
- class MDAnalysis.coordinates.base.ConverterBase[source]
Base class for converting to other libraries.
Deprecated since version 2.7.0: This class has been moved to
MDAnalysis.converters.base.ConverterBase
and will be removed fromMDAnalysis.coordinates.base
in 3.0.0.- 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.
- units = {'length': None, 'time': None, 'velocity': None}
dict with units of of time and length (and velocity, force, … for formats that support it)
7.33.5. Helper classes
The following classes contain basic functionality that all readers and writers share.
- class MDAnalysis.coordinates.base.IOBase[source]
Base class bundling common functionality for trajectory I/O.
Changed in version 0.8: Added context manager protocol.
- convert_forces_from_native(force, inplace=True)[source]
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)[source]
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)[source]
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)[source]
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)[source]
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)[source]
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)[source]
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)[source]
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.
- units = {'length': None, 'time': None, 'velocity': None}
dict with units of of time and length (and velocity, force, … for formats that support it)