6.23. FHI-AIMS file format — MDAnalysis.coordinates.FHIAIMS¶
Classes to read and write FHI-AIMS coordinate files.
The cell vectors are specified by the (optional) lines with the
lattice_vector tag:
lattice_vector x  y  z
where x, y, and z are expressed in ångström (Å).
Note
In the original FHI-AIMS format, up to three lines with
lattice_vector are allowed (order matters) where the absent line implies
no periodicity in that direction.  In MDAnalysis, only the case of no
lattice_vector or three lattice_vector lines are allowed.
Atomic positions and names are specified either by the atom or by the
atom_frac tags:
atom           x  y  z  name
atom_frac      nx ny nz name
where x, y, and z are expressed in ångström, and nx, ny and nz are real numbers in [0, 1] and are used to compute the atomic positions in units of the basic cell.
Atomic velocities can be added on the line right after the corresponding
atom in units of Å/ps using the velocity tag:
velocity      vx vy vz
The field name is a string identifying the atomic species. See also the specifications in the official FHI-AIMS format.
6.23.1. Classes¶
- 
class MDAnalysis.coordinates.FHIAIMS.Timestep(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 timeis 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.- 
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- Example - Using a Python - sliceobject:- 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. 
 - 
dimensions¶
- unitcell dimensions (A, B, C, alpha, beta, gamma) 
 - 
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. 
 - 
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. 
 - 
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. 
 - 
has_forces¶
- A boolean of whether this Timestep has force data - This can be changed to - Trueor- Falseto allocate space for or remove the data.- New in version 0.11.0. 
 - 
has_positions¶
- A boolean of whether this Timestep has position data - This can be changed to - Trueor- Falseto 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 - Trueor- Falseto allocate space for or remove the data.- New in version 0.11.0. 
 - 
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 
 - 
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 - NoDataErrorwhen no position data present
 - 
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. 
 - 
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. 
 - 
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. 
 - 
volume¶
- volume of the unitcell 
 
- 
class MDAnalysis.coordinates.FHIAIMS.FHIAIMSReader(filename, convert_units=True, n_atoms=None, **kwargs)[source]¶
- Reader for the FHIAIMS geometry format. - Single frame reader for the FHI-AIMS input file format. Reads geometry (3D and molecules only), positions (absolut or fractional), velocities if given, all according to the FHI-AIMS format specifications - 
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, n_atoms=None, **kwargs)[source]¶
- Returns a FHIAIMSWriter for filename. - Parameters: - filename (str) – filename of the output FHI-AIMS file - Returns: - Return type: - FHIAIMSWriter
 - 
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 - AuxReaderinstance, or the data itself as e.g. a filename; in the latter case an appropriate- AuxReaderis 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.auxnamespace 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.pullor- u.trajectory.ts.aux['pull'].- 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 - Timestepobject as argument, which will be transformed and returned to the Reader. The transformations can be part of the- transformationsmodule, 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 
 - 
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). Nonecorresponds to the default of 0, i.e., the initial frame.
- stop (int or None) – Last frame index (exclusive). Nonecorresponds 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, Nonecorresponds 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 a- slicewhen- stop=Noneand- 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.AnalysisBasethe values returned by- check_slice_indices()are used to splice the trajectory by creating a- sliceinstance- 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)] # [] 
- start (int or None) – Starting frame index (inclusive). 
 - 
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 - Falseso 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 - Falseso 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 - Falseso 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 - Falseso 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()¶
- 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: - 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: - 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. - 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.
- stop, step) ((start,) – Options for iterating over a slice of the auxiliary.
- selected (lst | ndarray, optional) – List of steps to iterate over.
 - Yields: - AuxStepobject- See also 
 - 
next()¶
- Forward one step to next frame. 
 - 
next_as_aux(auxname)¶
- Move to the next timestep for which there is at least one step from the auxiliary auxname within the cutoff specified in auxname. - This allows progression through the trajectory without encountering - NaNrepresentative 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 
 - 
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
 - 
rename_aux(auxname, new)¶
- Change the name of the auxiliary auxname to new. - Provided there is not already an auxiliary named new, the auxiliary name will be changed in ts.aux namespace, the trajectory’s list of added auxiliaries, and in the auxiliary reader itself. - Parameters: - Raises: - ValueError– If the name new is already in use by an existing auxiliary.
 - 
rewind()¶
- Position at beginning of trajectory 
 - 
set_aux_attribute(auxname, attrname, new)¶
- Set the value of attrname in the auxiliary auxname. - Parameters: - See also 
 - 
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.FHIAIMS.FHIAIMSWriter(filename, convert_units=True, n_atoms=None, **kwargs)[source]¶
- FHI-AIMS Writer. - Single frame writer for the FHI-AIMS format. Writes geometry (3D and molecules only), positions (absolut only), velocities if given, all according to the FHI-AIMS format specifications. - If no atom names are given, it will set each atom name to “X”. - Set up the FHI-AIMS Writer - Parameters: - 
close()¶
- Close the trajectory file. 
 - 
convert_dimensions_to_unitcell(ts, inplace=True)¶
- Read dimensions from timestep ts and return appropriate unitcell. - The default is to return - [A,B,C,alpha,beta,gamma]; if this is not appropriate then this method has to be overriden.
 - 
convert_forces_from_native(force, inplace=True)¶
- Conversion of forces array force from native to base units - Parameters: - force (array_like) – Forces to transform
- inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data
 - Note - By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided. - 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 - Falseso 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 - Falseso 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 - Falseso 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 - Falseso 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. 
 - 
fmt= {'box_triclinic': 'lattice_vector {box[0]:12.8f} {box[1]:12.8f} {box[2]:12.8f}\nlattice_vector {box[3]:12.8f} {box[4]:12.8f} {box[5]:12.8f}\nlattice_vector {box[6]:12.8f} {box[7]:12.8f} {box[8]:12.8f}\n', 'vel': 'velocity {vel[0]:12.8f} {vel[1]:12.8f} {vel[2]:12.8f}\n', 'xyz': 'atom {pos[0]:12.8f} {pos[1]:12.8f} {pos[2]:12.8f} {name:<3s}\n'}¶
- format strings for the FHI-AIMS file (all include newline) 
 - 
has_valid_coordinates(criteria, x)¶
- Returns - Trueif 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: 
 - 
write(obj)¶
- Write current timestep, using the supplied obj. - Parameters: - obj ( - AtomGroupor- 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. 
 
- 
6.23.2. Developer notes: FHIAIMSWriter format strings¶
The FHIAIMSWriter class has a FHIAIMSWriter.fmt
attribute, which is a dictionary of different strings for writing
lines in .in files.  These are as follows:
- xyz
- An atom line without velocities. Requires that the name and pos keys be supplied. E.g.: - fmt['xyz'].format(pos=(0.0, 1.0, 2.0), name='O') 
- vel
- An line that specifies velocities: - fmt['xyz'].format(vel=(0.1, 0.2, 0.3)) 
- box_triclinic
- The (optional) initial lines of the file which gives box dimensions. Requires the box keyword, as a length 9 vector. This is a flattened version of the (3, 3) triclinic vector representation of the unit cell.