AtomArrayStack
#
- class biotite.structure.AtomArrayStack(depth, length)[source]#
Bases:
_AtomArrayBase
A collection of multiple
AtomArray
instances, where each atom array has equal annotation arrays.Effectively, this means that each atom is occuring in every array in the stack at differing coordinates. This situation arises e.g. in NMR-elucidated or simulated structures. Since the annotations are equal for each array, the annotation arrays are 1-D, while the coordinate array is 3-D (m x n x 3). A detailed description of each annotation category can be viewed here.
Indexing works similar to
AtomArray
, with the difference, that two index dimensions are possible: The first index dimension specifies the array(s), the second index dimension specifies the atoms in each array (same as the index inAtomArray
). Using a single integer as first dimension index returns a singleAtomArray
instance.Concatenation of atoms for each array in the stack is done using the ‘+’ operator. For a list of
AtomArray
objects, useconcatenate()
. For addition of atom arrays onto the stack use thestack()
method.The
box
attribute has the shape m x 3 x 3, as the cell might be different for each frame in the atom array stack.- Parameters:
- depthint
The fixed amount of arrays in the stack. When indexing, this is the length of the first dimension.
- lengthint
The fixed amount of atoms in each array in the stack. When indexing, this is the length of the second dimension.
- Attributes:
- {annot}ndarray, shape=(n,)
Mutliple n-length annotation arrays.
- coordndarray, dtype=float, shape=(m,n,3)
ndarray containing the x, y and z coordinate of the atoms.
- bonds: BondList or None
A
BondList
, specifying the indices of atoms that form a chemical bond.- box: ndarray, dtype=float, shape=(m,3,3) or None
The surrounding box. May represent a MD simulation box or a crystallographic unit cell.
shape
tuple of intTuple of array dimensions.
See also
AtomArray
Representation of a single structure model.
Examples
Creating an atom array stack from two arrays:
>>> atom1 = Atom([1,2,3], chain_id="A") >>> atom2 = Atom([2,3,4], chain_id="A") >>> atom3 = Atom([3,4,5], chain_id="B") >>> atom_array1 = array([atom1, atom2, atom3]) >>> print(atom_array1.coord) [[1. 2. 3.] [2. 3. 4.] [3. 4. 5.]] >>> atom_array2 = atom_array1.copy() >>> atom_array2.coord += 3 >>> print(atom_array2.coord) [[4. 5. 6.] [5. 6. 7.] [6. 7. 8.]] >>> array_stack = stack([atom_array1, atom_array2]) >>> print(array_stack.coord) [[[1. 2. 3.] [2. 3. 4.] [3. 4. 5.]] [[4. 5. 6.] [5. 6. 7.] [6. 7. 8.]]]
- add_annotation(category, dtype)#
Add an annotation category, if not already existing.
Initially the new annotation is filled with the zero representation of the given type.
- Parameters:
- categorystr
The annotation category to be added.
- dtypetype or str
A type instance or a valid NumPy dtype string. Defines the type of the annotation.
See also
set_annotation
Assign directly a value to an annotation.
Notes
If the annotation category already exists, a compatible dtype is chosen, that is also able to represent the old values.
- array_length()#
Get the length of the atom array.
This value is equivalent to the length of each annotation array. For
AtomArray
it is the same aslen(array)
.- Returns:
- lengthint
Length of the array(s).
- copy()#
Create a deep copy of this object.
- Returns:
- copy
A copy of this object.
- del_annotation(category)#
Removes an annotation category.
- Parameters:
- categorystr
The annotation category to be removed.
- equal_annotation_categories(item)#
Check, if this object shares equal annotation array categories with the given
AtomArray
orAtomArrayStack
.- Parameters:
- itemAtomArray or AtomArrayStack
The object to compare the annotation arrays with.
- Returns:
- equalitybool
True, if the annotation array names are equal.
- equal_annotations(item, equal_nan=True)#
Check, if this object shares equal annotation arrays with the given
AtomArray
orAtomArrayStack
.- Parameters:
- itemAtomArray or AtomArrayStack
The object to compare the annotation arrays with.
- equal_nanbool
Whether to count nan values as equal. Default: True.
- Returns:
- equalitybool
True, if the annotation arrays are equal.
- get_annotation(category)#
Return an annotation array.
- Parameters:
- categorystr
The annotation category to be returned.
- Returns:
- arrayndarray
The annotation array.
- get_annotation_categories()#
Return a list containing all annotation array categories.
- Returns:
- categorieslist
The list containing the names of each annotation array.
- get_array(index)#
Obtain the atom array instance of the stack at the specified index.
The same as
stack[index]
, if index is an integer.- Parameters:
- indexint
Index of the atom array.
- Returns:
- arrayAtomArray
AtomArray at position index.
- set_annotation(category, array)#
Set an annotation array. If the annotation category does not exist yet, the category is created.
- Parameters:
- categorystr
The annotation category to be set.
- arrayndarray
The new value of the annotation category. The size of the array must be the same as the array length.
Notes
If the annotation category already exists, a compatible dtype is chosen, that is able to represent the old and new array values.
- stack_depth()#
Get the depth of the stack.
This value represents the amount of atom arrays in the stack. It is the same as
len(array)
.- Returns:
- lengthint
Length of the array(s).
Gallery#

Visualization of normal modes from an elastic network model