ITIM

Module: itim

class pytim.itim.ITIM(universe, group=None, alpha=1.5, normal='guess', molecular=True, max_layers=1, radii_dict=None, cluster_cut=None, cluster_threshold_density=None, extra_cluster_groups=None, info=False, centered=False, warnings=False, mesh=0.4, autoassign=True, **kargs)

Identifies interfacial molecules at macroscopically flat interfaces.

(Pártay, L. B.; Hantal, Gy.; Jedlovszky, P.; Vincze, Á.; Horvai, G., J. Comp. Chem. 29, 945, 2008)

Parameters:
  • universe (Object) – The MDAnalysis Universe, MDTraj trajectory or OpenMM Simulation objects.
  • group (Object) – An AtomGroup, or an array-like object with the indices of the atoms in the group. Will identify the interfacial molecules from this group
  • alpha (float) – The probe sphere radius
  • normal (str) – The macroscopic interface normal direction ‘x’,’y’, ‘z’ or ‘guess’ (default)
  • molecular (bool) – Switches between search of interfacial molecules / atoms (default: True)
  • max_layers (int) – The number of layers to be identified
  • radii_dict (dict) – Dictionary with the atomic radii of the elements in the group. If None is supplied, the default one (from GROMOS 43a1) will be used.
  • cluster_cut (float) – Cutoff used for neighbors or density-based cluster search (default: None disables the cluster analysis)
  • cluster_threshold_density (float) – Number density threshold for the density-based cluster search. ‘auto’ determines the threshold automatically. Default: None uses simple neighbors cluster search, if cluster_cut is not None
  • extra_cluster_groups (Object) – Additional groups, to allow for mixed interfaces
  • info (bool) – Print additional info
  • centered (bool) – Center the group
  • warnings (bool) – Print warnings
  • mesh (float) – The grid spacing used for the testlines (default 0.4 Angstrom)
  • autoassign (bool) – If true (default) detect the interface every time a new frame is selected.

Example:

>>> import MDAnalysis as mda
>>> import numpy as np
>>> import pytim
>>> from pytim.datafiles import *
>>>
>>> u         = mda.Universe(WATER_GRO)
>>> oxygens   = u.select_atoms("name OW")
>>>
>>> interface = pytim.ITIM(u, alpha=1.5, max_layers=4,molecular=True)
>>> # atoms in the layers can be accesses either through
>>> # the layers array:
>>> print (interface.layers)
[[<AtomGroup with 786 atoms> <AtomGroup with 681 atoms>
  <AtomGroup with 663 atoms> <AtomGroup with 651 atoms>]
 [<AtomGroup with 786 atoms> <AtomGroup with 702 atoms>
  <AtomGroup with 666 atoms> <AtomGroup with 636 atoms>]]
>>> interface.layers[0,0] # upper side, first layer
<AtomGroup with 786 atoms>
>>> interface.layers[1,2] # lower side, third layer
<AtomGroup with 666 atoms>
>>> # or as a whole AtomGroup. This can include all atoms in all layers
>>> interface.atoms
<AtomGroup with 5571 atoms>
>>> selection = interface.atoms.sides == 0
>>> interface.atoms[ selection ] # all atoms in the upper side layer
<AtomGroup with 2781 atoms>
>>> selection = np.logical_and(interface.atoms.layers == 2 , selection)
>>> interface.atoms[ selection ] # upper side, second layer
<AtomGroup with 681 atoms>
>>> # the whole system can be quickly saved to a pdb file
>>> # including the layer information, written in the beta field
>>> # using:
>>> interface.writepdb('system.pdb',centered=True)
>>> # of course, the native interface of MDAnalysis can be used to
>>> # write pdb files, but the centering options are not available.
>>> # Writing to other formats that do not support the beta factor
>>> # will loose the information on the layers.
>>> interface.atoms.write('only_layers.pdb')
>>> # In some cases it might be necessary to compute two interfaces.
>>> # This could be done in the following way:
>>> import MDAnalysis as mda
>>> import pytim
>>> from pytim.datafiles import WATER_GRO, WATER_XTC
>>> u = mda.Universe(WATER_GRO,WATER_XTC)
>>> u2 = mda.Universe(WATER_GRO,WATER_XTC)
>>> inter = pytim.ITIM(u,group=u.select_atoms('resname SOL'))
>>> inter2 = pytim.ITIM(u2,group=u2.select_atoms('resname SOL'))
>>> for ts in u.trajectory[::50]:
...     ts2 = u2.trajectory[ts.frame]
autoassign

(bool) assign layers every time a frame changes

label_planar_sides()

Assign to all layers a label (the beta tempfactor) that can be used in pdb files. Additionally, set the new layers and sides.

layers

Access the layers as numpy arrays of AtomGroups.

The object can be sliced as usual with numpy arrays, so, for example:

>>> import MDAnalysis as mda
>>> import pytim
>>> from pytim.datafiles import *
>>>
>>> u         = mda.Universe(WATER_GRO)
>>> oxygens   = u.select_atoms("name OW")
>>>
>>> interface = pytim.ITIM(u, alpha=1.5, max_layers=4,molecular=True)
>>> print(interface.layers[0,:])  # upper side (0), all layers
[<AtomGroup with 786 atoms> <AtomGroup with 681 atoms>
 <AtomGroup with 663 atoms> <AtomGroup with 651 atoms>]
>>> repr(interface.layers[1,0])  # lower side (1), first layer (0)
'<AtomGroup with 786 atoms>'
>>> print(interface.layers[:,0:3]) # 1st - 3rd layer (0:3), on both sides
[[<AtomGroup with 786 atoms> <AtomGroup with 681 atoms>
  <AtomGroup with 663 atoms>]
 [<AtomGroup with 786 atoms> <AtomGroup with 702 atoms>
  <AtomGroup with 666 atoms>]]
>>> print(interface.layers[1,0:4:2]) # side 1, layers 1-4 & stride 2 (0:4:2)
[<AtomGroup with 786 atoms> <AtomGroup with 666 atoms>]
writepdb(filename='layers.pdb', centered='no', group='all', multiframe=True, tempfactors=None)

Write the frame to a pdb file, marking the atoms belonging to the layers with different beta factors.

Parameters:
  • filename (str) – the output file name
  • centered (str) – ‘origin’, ‘middle’, or ‘no’
  • group (AtomGroup) – if ‘all’ is passed, use universe
  • multiframe (bool) – append to pdb file if True
  • tempfactors (ndarray) – use this array as temp (beta) factors
Example: save the positions (centering the interface in the cell)
without appending
>>> import pytim
>>> import pytim.datafiles
>>> import MDAnalysis as mda
>>> from pytim.datafiles import WATER_GRO
>>> u = mda.Universe(WATER_GRO)
>>> interface = pytim.ITIM(u)
>>> interface.writepdb('layers.pdb',multiframe=False)
Example: save the positions without centering the interface. This
will not shift the atoms from the original position (still, they will be put into the basic cell). The multiframe option set to False will overwrite the file
>>> interface.writepdb('layers.pdb',centered='no')

Note that if GITIM is used, and the symmetry option is different from 'planar', the centered='origin' option is equivalent to centered='middle'.