6.27. Base classes — MDAnalysis.coordinates.base

Derive other Timestep, FrameIterator, Reader and Writer classes from the classes in this module. The derived classes must follow the Trajectory API in MDAnalysis.coordinates.__init__.

6.27.1. Timestep

A Timestep holds information for the current time frame in the trajectory. It is one of the central data structures in MDAnalysis.

class MDAnalysis.coordinates.base.Timestep[source]
__init__(n_atoms, **kwargs)[source]

Create a Timestep, representing a frame of a trajectory

Parameters:
  • n_atoms (int) – The total number of atoms this Timestep describes
  • positions (bool, optional) – Whether this Timestep has position information [True]
  • velocities (bool (optional)) – Whether this Timestep has velocity information [False]
  • forces (bool (optional)) – Whether this Timestep has force information [False]
  • reader (Reader (optional)) – A weak reference to the owning Reader. Used for when attributes require trajectory manipulation (e.g. dt)
  • dt (float (optional)) – The time difference between frames (ps). If time is set, then dt will be ignored.
  • time_offset (float (optional)) – The starting time from which to calculate time (in ps)

Changed in version 0.11.0: Added keywords for positions, velocities and forces. Can add and remove position/velocity/force information by using the has_* attribute.

classmethod from_coordinates(positions=None, velocities=None, forces=None, **kwargs)[source]

Create an instance of this Timestep, from coordinate data

Can pass position, velocity and force data to form a Timestep.

New in version 0.11.0.

classmethod from_timestep(other, **kwargs)[source]

Create a copy of another Timestep, in the format of this Timestep

New in version 0.11.0.

_init_unitcell()[source]

Create custom datastructure for _unitcell.

n_atoms

A read only view of the number of atoms this Timestep has

Changed in version 0.11.0: Changed to read only property

time

The time in ps of this timestep

This is calculated as:

time = ts.data['time_offset'] + ts.time

Or, if the trajectory doesn’t provide time information:

time = ts.data['time_offset'] + ts.frame * ts.dt

New in version 0.11.0.

dt

The time difference in ps between timesteps

Note

This defaults to 1.0 ps in the absence of time data

New in version 0.11.0.

positions

A record of the positions of all atoms in this Timestep

Setting this attribute will add positions to the Timestep if they weren’t originally present.

Returns:positions – position data of shape (n_atoms, 3) for all atoms
Return type:numpy.ndarray with dtype numpy.float32
Raises:MDAnalysis.exceptions.NoDataError – if the Timestep has no position data

Changed in version 0.11.0: Now can raise NoDataError when no position data present

velocities

A record of the velocities of all atoms in this Timestep

Setting this attribute will add velocities to the Timestep if they weren’t originally present.

Returns:velocities – velocity data of shape (n_atoms, 3) for all atoms
Return type:numpy.ndarray with dtype numpy.float32
Raises:MDAnalysis.exceptions.NoDataError – if the Timestep has no velocity data

New in version 0.11.0.

forces

A record of the forces of all atoms in this Timestep

Setting this attribute will add forces to the Timestep if they weren’t originally present.

Returns:forces – force data of shape (n_atoms, 3) for all atoms
Return type:numpy.ndarray with dtype numpy.float32
Raises:MDAnalysis.exceptions.NoDataError – if the Timestep has no force data

New in version 0.11.0.

has_positions

A boolean of whether this Timestep has position data

This can be changed to True or False to allocate space for or remove the data.

New in version 0.11.0.

has_velocities

A boolean of whether this Timestep has velocity data

This can be changed to True or False to allocate space for or remove the data.

New in version 0.11.0.

has_forces

A boolean of whether this Timestep has force data

This can be changed to True or False to allocate space for or remove the data.

New in version 0.11.0.

_pos

numpy.ndarray of dtype float32 of shape (n_atoms, 3) and internal FORTRAN order, holding the raw cartesian coordinates (in MDAnalysis units, i.e. Å).

Note

Normally one does not directly access _pos but uses the coordinates() method of an AtomGroup but sometimes it can be faster to directly use the raw coordinates. Any changes to this array are immediately reflected in atom positions. If the frame is written to a new trajectory then the coordinates are changed. If a new trajectory frame is loaded, then all contents of _pos are overwritten.

_velocities

numpy.ndarray of dtype float32. of shape (n_atoms, 3), holding the raw velocities (in MDAnalysis units, i.e. typically Å/ps).

Note

Normally velocities are accessed through the velocities or the velocities() method of an AtomGroup

_velocities only exists if the has_velocities flag is True

New in version 0.7.5.

_forces

numpy.ndarray of dtype float32. of shape (n_atoms, 3), holding the forces

_forces only exists if has_forces is True

New in version 0.11.0: Added as optional to Timestep

dimensions

unitcell dimensions (A, B, C, alpha, beta, gamma)

lengths a, b, c are in the MDAnalysis length unit (Å), and angles are in degrees.

Setting dimensions will populate the underlying native format description of box dimensions

triclinic_dimensions

The unitcell dimensions represented as triclinic vectors

Returns:A (3, 3) numpy.ndarray of unit cell vectors
Return type:numpy.ndarray

Examples

The unitcell for a given system can be queried as either three vectors lengths followed by their respective angle, or as three triclinic vectors.

>>> ts.dimensions
array([ 13.,  14.,  15.,  90.,  90.,  90.], dtype=float32)
>>> ts.triclinic_dimensions
array([[ 13.,   0.,   0.],
       [  0.,  14.,   0.],
       [  0.,   0.,  15.]], dtype=float32)

Setting the attribute also works:

>>> ts.triclinic_dimensions = [[15, 0, 0], [5, 15, 0], [5, 5, 15]]
>>> ts.dimensions
array([ 15.        ,  15.81138802,  16.58312416,  67.58049774,
        72.45159912,  71.56504822], dtype=float32)

New in version 0.11.0.

volume

volume of the unitcell

__getitem__(atoms)[source]

Get a selection of coordinates

ts[i]

return coordinates for the i’th atom (0-based)

ts[start:stop:skip]

return an array of coordinates, where start, stop and skip correspond to atom indices, MDAnalysis.core.groups.Atom.index (0-based)
__eq__(other)[source]

Compare with another Timestep

New in version 0.11.0.

__iter__()[source]

Iterate over coordinates

for x in ts

iterate of the coordinates, atom by atom
copy()[source]

Make an independent (“deep”) copy of the whole Timestep.

copy_slice(sel)[source]

Make a new Timestep containing a subset of the original Timestep.

Parameters:sel (array_like or slice) – The underlying position, velocity, and force arrays are sliced using a list, slice, or any array-like.
Returns:A Timestep object of the same type containing all header information and all atom information relevant to the selection.
Return type:Timestep

Note

The selection must be a 0 based slice or array of the atom indices in this Timestep

Example

Using a Python slice object:

new_ts = ts.copy_slice(slice(start, stop, step))

Using a list of indices:

new_ts = ts.copy_slice([0, 2, 10, 20, 23])

New in version 0.8.

Changed in version 0.11.0: Reworked to follow new Timestep API. Now will strictly only copy official attributes of the Timestep.

6.27.2. FrameIterators

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.

New 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

FrameIteratorBase,

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

FrameIteratorBase,

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.

6.27.3. 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:

  1. Readers for multi frame trajectories, i.e., file formats that typically contain many frames. These readers are typically derived from ReaderBase.
  2. Readers for single frame formats: These file formats only contain a single coordinate set. These readers are derived from :class`: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 a close() method, to ensure proper release of resources (mainly file handles). Readers that are inherently safe in this regard should subclass ProtoReader instead.

See the Trajectory API definition in MDAnalysis.coordinates.__init__ for the required attributes and methods.

See also

ProtoReader

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 now super() through this class. Added attribute _ts_kwargs, which is created in init. Provides kwargs to be passed to Timestep

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(auxname, auxdata, format=None, **kwargs)

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 appropriate AuxReader 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 with trajectory[i]).

The representative value(s) of the auxiliary data for each timestep (as calculated by the AuxReader) are stored in the current timestep in the ts.aux namespace under auxname; 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 or u.trajectory.ts.aux['pull'].

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 the transformations 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 to

for 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
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:

tuple (int, int, int)

Warning

The returned values start, stop and step give the expected result when passed in range() but gives unexpected behavior when passed in a slice when stop=None and step=-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() to range()

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 by check_slice_indices() are used to splice the trajectory by creating a slice instance

slice(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.

New 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.

New 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.

New 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.

New 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

dt

Time between two trajectory frames in picoseconds.

frame

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:
  • auxname (str) – Name of the auxiliary to get value for
  • attrname (str) – Name of gettable attribute in the auxiliary reader
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:list
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.

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.
  • stop, step) ((start,) – Options for iterating over a slice of the auxiliary.
  • selected (lst | ndarray, optional) – List of steps to iterate over.
Yields:

AuxStep object

See also

iter_as_aux()

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 calling next().

See the Auxiliary API.

See also

iter_as_aux()

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:int
Raises:NotImplementedError – when the number of atoms can’t be deduced
remove_auxiliary(auxname)

Clear data and close the AuxReader for the auxiliary auxname.

See also

add_auxiliary()

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:
  • auxname (str) – Name of the auxiliary to rename
  • new (str) – New name to try set
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:
  • auxname (str) – Name of the auxiliary to alter
  • attrname (str) – Name of settable attribute in the auxiliary reader
  • new – New value to try set attrname to
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

totaltime

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.

transformations

Returns the list of transformations

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 a Timestep object.

New 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

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(auxname, auxdata, format=None, **kwargs)

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 appropriate AuxReader 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 with trajectory[i]).

The representative value(s) of the auxiliary data for each timestep (as calculated by the AuxReader) are stored in the current timestep in the ts.aux namespace under auxname; 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 or u.trajectory.ts.aux['pull'].

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 the transformations 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
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:

tuple (int, int, int)

Warning

The returned values start, stop and step give the expected result when passed in range() but gives unexpected behavior when passed in a slice when stop=None and step=-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() to range()

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 by check_slice_indices() are used to splice the trajectory by creating a slice instance

slice(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()[source]

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.

New 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.

New 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.

New 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.

New 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

dt

Time between two trajectory frames in picoseconds.

frame

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:
  • auxname (str) – Name of the auxiliary to get value for
  • attrname (str) – Name of gettable attribute in the auxiliary reader
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:list
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.

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.
  • stop, step) ((start,) – Options for iterating over a slice of the auxiliary.
  • selected (lst | ndarray, optional) – List of steps to iterate over.
Yields:

AuxStep object

See also

iter_as_aux()

next()[source]

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 calling next().

See the Auxiliary API.

See also

iter_as_aux()

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:int
Raises:NotImplementedError – when the number of atoms can’t be deduced
remove_auxiliary(auxname)

Clear data and close the AuxReader for the auxiliary auxname.

See also

add_auxiliary()

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:
  • auxname (str) – Name of the auxiliary to rename
  • new (str) – New name to try set
Raises:

ValueError – If the name new is already in use by an existing auxiliary.

rewind()[source]

Position at beginning of trajectory

set_aux_attribute(auxname, attrname, new)

Set the value of attrname in the auxiliary auxname.

Parameters:
  • auxname (str) – Name of the auxiliary to alter
  • attrname (str) – Name of settable attribute in the auxiliary reader
  • new – New value to try set attrname to
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

totaltime

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.

transformations

Returns the list of transformations

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

ReaderBase

Changed in version 0.11.0: Frames now 0-based instead of 1-based

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(auxname, auxdata, format=None, **kwargs)[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 appropriate AuxReader 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 with trajectory[i]).

The representative value(s) of the auxiliary data for each timestep (as calculated by the AuxReader) are stored in the current timestep in the ts.aux namespace under auxname; 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 or u.trajectory.ts.aux['pull'].

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 the transformations 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 to

for 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
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:

tuple (int, int, int)

Warning

The returned values start, stop and step give the expected result when passed in range() but gives unexpected behavior when passed in a slice when stop=None and step=-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() to range()

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 by check_slice_indices() are used to splice the trajectory by creating a slice instance

slice(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)]  # []
dt

Time between two trajectory frames in picoseconds.

frame

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:
  • auxname (str) – Name of the auxiliary to get value for
  • attrname (str) – Name of gettable attribute in the auxiliary reader
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:list
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.

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.
  • stop, step) ((start,) – Options for iterating over a slice of the auxiliary.
  • selected (lst | ndarray, optional) – List of steps to iterate over.
Yields:

AuxStep object

See also

iter_as_aux()

next()[source]

Forward one step to next frame.

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 calling next().

See the Auxiliary API.

See also

iter_as_aux()

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:int
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

add_auxiliary()

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:
  • auxname (str) – Name of the auxiliary to rename
  • new (str) – New name to try set
Raises:

ValueError – If the name new is already in use by an existing auxiliary.

rewind()[source]

Position at beginning of trajectory

set_aux_attribute(auxname, attrname, new)[source]

Set the value of attrname in the auxiliary auxname.

Parameters:
  • auxname (str) – Name of the auxiliary to alter
  • attrname (str) – Name of settable attribute in the auxiliary reader
  • new – New value to try set attrname to
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

totaltime

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.

transformations

Returns the list of transformations

6.27.4. 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 MDAnalysis.coordinates.__init__ for the required attributes and methods.

Deprecated since version 1.0.0: write_next_timestep() has been deprecated, please use write() 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.

New 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.

New 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.

New 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.

New 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
Returns:

Return type:

bool

write(obj)[source]

Write current timestep, using the supplied obj.

Parameters:obj (AtomGroup or Universe) – write coordinate information associate with obj

Note

The size of the obj must be the same as the number of atoms provided when setting up the trajectory.

Deprecated since version 1.0.0: Deprecated the use of Timestep as arguments to write. Use either an AtomGroup or Universe. To be removed in version 2.0.

write_next_timestep(obj)[source]

Write current timestep, using the supplied obj.

Parameters:obj (AtomGroup or Universe) – write coordinate information associated with obj

Deprecated since version 1.0.0: Deprecated, use write() instead

6.27.5. Converters

Converters output information to other libraries.

class MDAnalysis.coordinates.base.ConverterBase[source]

Base class for converting to other libraries.

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.

New 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.

New 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.

New 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.

New in version 0.7.5.

6.27.6. 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.

close()[source]

Close the trajectory file.

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.

New 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.

New 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.

New 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.

New 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)