SimpleInterface

class pytim.simple_interface.SimpleInterface(universe, group=None, alpha=1.5, symmetry='generic', normal='z', include_zero_radius=False, upper=None, lower=None)[source]

This simple interface is designed to allow the user to define its own interface by supplying one or more groups, and be able to use all other analysis tools of pytim, which require an interface, also when it is not advisable/possible to use any of the other classes

Parameters:
  • universe (Universe) – The MDAnalysis Universe.

  • group (AtomGroup) – The AtomGroup with the interfacial atoms. (only for symmetry=’generic’). See the upper and lower groups options below.

  • alpha (float) – The probe sphere radius (used for handling the periodicity in the surface layer triangulation only)

  • normal (str) – The macroscopic interface normal direction ‘x’,’y’, or ‘z’(default)

  • symmetry (str) – Either ‘grneric’ (default) or ‘planar’. Selects the type of interpolation/distance calculation. If ‘planar’ two groups have to be also passed (upper and lower)

  • upper (AtomGroup) – The upper surface group (if symmetry=’planar’)

  • lower (AtomGroup) – The lower surface group (if symmetry=’planar’)

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

Example: compute an intrinsic profile by interpolating the surface position from the P atoms of a DPPC bilayer

>>> import pytim
>>> import numpy as np
>>> from pytim.datafiles import DPPC_GRO
>>> import MDAnalysis as mda
>>> u = mda.Universe(DPPC_GRO)
>>> p = u.select_atoms('name P8')
>>>
>>> dppc  = u.select_atoms('resname DPPC')
>>> water = u.select_atoms('resname SOL')
>>>
>>> box = u.dimensions[:3]
>>> # here we want to use as 'surface' atoms only the P atoms
>>> # and compute an intrinsic profile with respect to the
>>> # interpolated surface.
>>> # We need to distinguish the upper and lower leaflet by hand!
>>> upper = p[p.positions[:,2]>box[2]/2.]
>>> lower = p - upper
>>>
>>> # if symmetry=='planar', supply an upper and a lower group
>>> # if symmetry=='generic', supply only group
>>>
>>> # here we need a large value of alpha (for pbc reconstruction
>>> # of the triangulation used for the interpolation)
>>> inter = pytim.SimpleInterface(u,symmetry='planar', upper=upper, lower=lower,alpha=5.)
>>>
>>> profile1 = pytim.observables.Profile(interface=inter)
>>> profile2 = pytim.observables.Profile(interface=inter)
>>> np.random.seed(1)
>>> profile1.sample(dppc)
>>> profile2.sample(water)
>>>
>>> lo,up,av1 = profile1.get_values(binwidth=.5)
>>> lo,up,av2 = profile2.get_values(binwidth=.5)
>>> np.set_printoptions(8)
>>> print (av2[64:70])
[0.04282182 0.04154116 0.04147182        inf 0.04651406 0.05234671]
property include_zero_radius

(bool) include atoms with zero radius in the analysis (excluded by default)

property analysis_group

(AtomGroup) the group, the surface of which should be computed

property autoassign

(bool) assign layers every time a frame changes

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.

property layers

AtomGroups of atoms in layers

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

property surface_cluster_cut

(float) distance cut to calculate the surface clusters

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