GITIM¶
Module: gitim¶
- class pytim.gitim.GITIM(universe, group=None, alpha=2.0, normal='guess', molecular=True, max_layers=1, radii_dict=None, cluster_cut=None, include_zero_radius=False, cluster_threshold_density=None, extra_cluster_groups=None, biggest_cluster_only=False, surface_cluster_cut=None, symmetry='generic', centered=False, info=False, warnings=False, autoassign=True, _noextrapoints=False, **kargs)[source]¶
Identifies interfacial molecules at curved interfaces.
(Sega, M.; Kantorovich, S.; Jedlovszky, P.; Jorge, M., J. Chem. Phys. 138, 044110, 2013)
- 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) – ‘x’,’y,’z’ or ‘guess’ (for planar interfaces only)
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
biggest_cluster_only (bool) – Tag as surface atoms/molecules only those in the largest cluster. Need to specify also a
cluster_cut
value.surface_cluster_cut (float) – Filter surface atoms/molecules to include only those in the largest cluster of (initially detected) surface ones. (default: None disables the filtering)
symmetry (str) – Gives the code a hint about the topology
symmetry – Gives the code a hint about the topology of the interface: ‘generic’ (default) or ‘planar’
centered (bool) – Center the
group
include_zero_radius (bool) – if false (default) exclude atoms with zero radius from the surface analysis (they are always included in the cluster search, if present in the relevant group) to avoid some artefacts.
info (bool) – Print additional info
warnings (bool) – Print warnings
autoassign (bool) – If true (default) detect the interface every time a new frame is selected.
Example:
>>> import MDAnalysis as mda >>> import pytim >>> from pytim.datafiles import * >>> >>> u = mda.Universe(MICELLE_PDB) >>> g = u.select_atoms('resname DPC') >>> >>> interface =pytim.GITIM(u,group=g,molecular=False, alpha=2.0) >>> layer = interface.layers[0] >>> interface.writepdb('gitim.pdb',centered=False) >>> print (repr(layer)) <AtomGroup with 909 atoms>
Successive layers can be identified with
GITIM
as well. In this example we identify two solvation shells of glucose:>>> import MDAnalysis as mda >>> import pytim >>> from pytim.datafiles import * >>> >>> u = mda.Universe(GLUCOSE_PDB) >>> g = u.select_atoms('name OW') >>> # it is faster to consider only oxygens. >>> # Hydrogen atoms are anyway within Oxygen's radius, >>> # in SPC* models. >>> interface =pytim.GITIM(u, group=g, alpha=2.0, max_layers=2) >>> >>> interface.writepdb('glucose_shells.pdb') >>> print (repr(interface.layers[0])) <AtomGroup with 54 atoms> >>> print (repr(interface.layers[1])) <AtomGroup with 117 atoms>
- property autoassign¶
(bool) assign layers every time a frame changes
- property surface_cluster_cut¶
(float) distance cut to calculate the surface clusters
- property include_zero_radius¶
(bool) include atoms with zero radius in the analysis (excluded by default)
- property layers¶
Access the layers as numpy arrays of AtomGroups.
The object can be sliced as usual with numpy arrays. Differently from
ITIM
, there are no sides. Example:>>> import MDAnalysis as mda >>> import pytim >>> from pytim.datafiles import MICELLE_PDB >>> >>> u = mda.Universe(MICELLE_PDB) >>> micelle = u.select_atoms('resname DPC') >>> inter = pytim.GITIM(u, group=micelle, max_layers=3,molecular=False) >>> inter.layers #all layers array([<AtomGroup with 909 atoms>, <AtomGroup with 301 atoms>, <AtomGroup with 164 atoms>], dtype=object) >>> inter.layers[0] # first layer (0) <AtomGroup with 909 atoms>
- property analysis_group¶
(AtomGroup) the group, the surface of which should be computed
- is_buried(pos)¶
Checks wether an array of positions are located below the first interfacial layer
- 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.
- prepare_box()¶
Before the analysis, pack every molecule into the box. Keep the original positions for latter use.
- reset_labels()¶
Reset labels before interfacial analysis
- 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 toFalse
will overwrite the file
>>> interface.writepdb('layers.pdb',centered='no')
Note that if
GITIM
is used, and thesymmetry
option is different fromplanar
, thecentered='origin'
option is equivalent tocentered='middle'
.