4.4.1. Generation and Analysis of HOLE pore profiles (Deprecated) — MDAnalysis.analysis.hole

Author:Lukas Stelzl, Oliver Beckstein
Year:2011-2012
Copyright:GNU Public License v2

With the help of this module, the hole program from the HOLE suite of tools [Smart1993] [Smart1996] can be run on frames in an MD trajectory or NMR ensemble in order to analyze an ion channel pore or transporter pathway [Stelzl2014] as a function of time or arbitrary order parameters. Data can be combined and analyzed.

HOLE must be installed separately and can be obtained in binary form from http://www.holeprogram.org/ or as source from https://github.com/osmart/hole2. (HOLE is open source and available under the Apache v2.0 license.)

4.4.1.1. Examples for using HOLE

The two classes HOLE and HOLEtraj in this module primarily act as wrappers around the hole program from the HOLE suite of tools. They contain many options to set options of hole. However, the defaults often work well and the following examples should be good starting points for applying HOLE to other problems.

4.4.1.1.1. Single structure

The following example runs hole on the experimental structure of the Gramicidin A (gA) channel. We use the HOLE class, which acts as a wrapper around hole. It therefore shares some of the limitations of HOLE, namely, that it can only process PDB files [1].

from MDAnalysis.analysis.hole import HOLE
from MDAnalysis.tests.datafiles import PDB_HOLE

H = HOLE(PDB_HOLE, executable="~/hole2/exe/hole")  # set path to your hole binary
H.run()
H.collect()
H.plot(linewidth=3, color="black", label=False)

The example assumes that hole was installed as ~/hole2/exe/hole)

4.4.1.1.2. Trajectory

One can also run hole on frames in a trajectory with HOLEtraj. In this case, provide a Universe:

import MDAnalysis as mda
from MDAnalysis.analysis.hole import HOLEtraj
from MDAnalysis.tests.datafiles import MULTIPDB_HOLE

u = mda.Universe(MULTIPDB_HOLE)
H = HOLEtraj(u, executable="~/hole2/exe/hole")
H.run()
H.plot3D()

The profiles are available as the attribute HOLEtraj.profiles (H.profiles in the example) and are indexed by frame number but can also be indexed by an arbitrary order parameter as shown in the next example.

4.4.1.1.3. Trajectory with RMSD as order parameter

In order to classify the HOLE profiles \(R(\zeta)\) the RMSD \(\rho\) to a reference structure is calculated for each trajectory frame (e.g. using the MDAnalysis.analysis.rms.RMSD analysis class). Then the HOLE profiles \(R_\rho(\zeta)\) can be ordered by the RMSD, which acts as an order parameter \(\rho\).

import MDAnalysis as mda
from MDAnalysis.analysis.hole import HOLEtraj
from MDAnalysis.analysis.rms import RMSD

from MDAnalysis.tests.datafiles import PDB_HOLE, MULTIPDB_HOLE

mda.start_logging()
ref = mda.Universe(PDB_HOLE)    # reference structure
u = mda.Universe(MULTIPDB_HOLE) # trajectory

# calculate RMSD
R = RMSD(u, reference=ref, select="protein", weights='mass')
R.run()

# HOLE analysis with order parameters
H = HOLEtraj(u, orderparameters=R.rmsd[:,2],
             executable="~/hole2/exe/hole")
H.run()

The HOLEtraj.profiles dictionary will have the order parameter as key for each frame. The plot functions will automatically sort the profiles by ascending order parameter. To access the individual profiles one can simply iterate over the sorted profiles (see HOLEtraj.sorted_profiles_iter()). For example, in order to plot the minimum radius as function of order parameter

\[r(\rho) = \min_\zeta R_\rho(\zeta)\]

we iterate over the profiles and process them in turn:

import numpy as np
import matplotlib.pyplot as plt

r_rho = np.array([[rho, profile.radius.min()] for rho, profile in H])

ax = plt.subplot(111)
ax.plot(r_rho[:, 0], r_rho[:, 1], lw=2)
ax.set_xlabel(r"order parameter RMSD $\rho$ ($\AA$)")
ax.set_ylabel(r"minimum HOLE pore radius $r$ ($\AA$)")

(The graph shows that the pore opens up with increasing order parameter.)

4.4.1.2. Data structures

A profile is stored as a numpy.recarray:

frame   rxncoord   radius
frame
integer frame number (only important when HOLE itself reads a trajectory)
rxncoord
the distance along the pore axis, in Å
radius
the pore radius, in Å

The HOLE.profiles or HOLEtraj.profiles dictionary holds one profile for each key. By default the keys are the frame numbers but HOLEtraj can take the optional orderparameters keyword argument and load an arbitrary order parameter for each frame. In that case, the key becomes the orderparameter.

Notes

The profiles dict is not ordered and hence one typically needs to manually order the keys first. Furthermore, duplicate keys are not possible: In the case of duplicate orderparameters, the last one read will be stored in the dict.

4.4.1.3. Analysis

class MDAnalysis.analysis.hole.HOLE(filename, **kwargs)[source]

Run hole on a single frame or a DCD trajectory.

hole is part of the HOLE suite of programs. It is used to analyze channels and cavities in proteins, especially ion channels.

Only a subset of all HOLE control parameters is supported and can be set with keyword arguments.

hole (as a FORTRAN77 program) has a number of limitations when it comes to filename lengths (must be shorter than the empirically found HOLE.HOLE_MAX_LENGTH). This class tries to work around them by creating temporary symlinks to files when needed but this can still fail when permissions are not correctly set on the current directory.

Running hole with the HOLE class is a 3-step process:

  1. set up the class with all desired parameters
  2. run hole with HOLE.run()
  3. collect the data from the output file with HOLE.collect()

The class also provides some simple plotting functions of the collected data such as HOLE.plot() or HOLE.plot3D().

New in version 0.7.7.

Changed in version 0.16.0: Added raseed keyword argument.

Set up parameters to run HOLE on PDB filename.

Parameters:
  • filename (string) – The filename is used as input for HOLE in the “COORD” card of the input file. It specifies the name of a PDB co-ordinate file to be used. This must be in Brookhaven protein databank format or something closely approximating this. Both ATOM and HETATM records are read. Note that if water molecules or ions are present in the channel these can be ignored on read by the use of the ignore_residues keyword. Wildcard pattern. A new feature (in release 2.1 of HOLE) was the option to include a wild card (*) in the filename. e.g., filename="ab*.pdb" will apply hole to all files in the directory whose name starts with ab and ends with .pdb. This is intended to aid the analysis of multiple copies of the same molecule - produced during molecular dynamics or other method. The hole procedure will be applied to each file in turn with the same setup conditions (initial point, sampling distance etc.). Graphics files will contain a combination of the individual runs, one after another. Note that the pdb files are read independently so that they need not have an identical number of atoms or atom order etc. (though they should be sufficiently similar for a HOLE run from identical starting conditions to be useful).
  • dcd (string, optional) – File name of DCD trajectory (must be supplied together with a matching PDB file filename) and then HOLE runs its analysis on each frame. It does multiple HOLE runs on positions taken from a CHARMM binary dynamics format DCD trajectory file. The dcd file must have exactly the same number of atoms in exactly the same order as the pdb file specified by filename. Note that if this option is used the pdb file is used as a template only - the coordinates are ignored. Note that structural parameters determined for each individual structure are written in a tagged format so that it is possible to extract the information from the text output file using a grep command. The reading of the file can be controlled by the step keyword and/or setting HOLE.dcd_iniskip to the number of frames to be skipped initially.
  • logfile (string, optional) – file name of the file collecting HOLE’s output (which can be parsed using HOLE.collect()); the default is “hole.out”.
  • sphpdb (string, optional) – path to the HOLE sph file, a PDB-like file containig the coordinates of the pore centers; the default is “hole.sph”. This keyword specifies the filename for output of the sphere centre information in pdb form. Its typical suffix is “.sph”. The co-ordinates are set to the sphere centres and the occupancies are the sphere radii. All centres are assigned the atom name QSS and residue name SPH and the residue number is set to the storage number of the centre. The file can be imported into molecular graphics programs but are likely to be bonded together in a awful manner - as the points are very close to one another. In VMD sph objects are best displayed as “Points”. Displaying .sph objects rather than rendered or dot surfaces can be useful to analyze the distance of particular atoms from the sphere-centre line. Most usefully .sph files can be used to produce molecular graphical output from a hole run. This is achieved by using the sph_process program to read the .sph file.
  • step (int, optional) – step size for going through the trajectory (skips step - 1 frames); the default is 1.
  • cpoint (array_like, optional) – coordinates of a point inside the pore, e.g. [12.3, 0.7, 18.55]. If set to None (the default) then HOLE’s own search algorithm is used. cpoint specifies a point which lies within the channel. For simple channels such as gramicidin results do not show great sensitivity to the exact point taken. An easy way to produce an initial point is to use molecular graphics to find two atoms which lie either side of the pore and to average their co-ordinates. Or if the channel structure contains water molecules or counter ions then take the coordinates of one of these (and use the ignore_residues keyword to ignore them in the pore radius calculation). If this card is not specified then HOLE now (from version 2.2) attempts to make a guess where the channel will be. The procedure assumes the channel is reasonably symmetric. The initial guess on cpoint will be the centroid of all alpha carbon atoms (name ‘CA’ in pdb file). This is then refined by a crude grid search up to 5 Å from the original position. This procedure works most of the time but is clearly far from infallible — results should be careful checked (with molecular graphics) if it is used.
  • cvect (array_like, optional) – Search direction, should be parallel to the pore axis, e.g. [0,0,1] for the z-axis. The default is None so that HOLE’s built-in procedure is used. If this keyword is None then HOLE attempts to make a guess where the channel will be. The procedure assumes the channel is reasonably symmetric. The guess will be either along the X axis (1,0,0), Y axis (0,1,0) or Z axis (0,0,1). If the structure is not aligned on one of these axis the results will clearly be approximate. If a guess is used then results should be carefully checked.
  • sample (float, optional) – distance of sample points in Å; the default is 0.2 Å. Specifies the distance between the planes used in the HOLE procedure. The default value should be reasonable for most purposes. However, if you wish to visualize a very tight constriction then specify a smaller value. This value determines how many points in the pore profile are calculated.
  • dotden (int, optional) – density of facettes for generating a 3D pore representation; default is 15. This number controls the density of dots which will be used by the program. A sphere of dots is placed on each centre determined in the Monte Carlo procedure. Only dots which do not lie within any other sphere are considered. The actual number of dots written is therefore controlled by dotden and sample. dotden should be set to between 5 (few dots per sphere) and 35 (large number of dots per sphere).
  • endrad (float, optional) – Radius which is considered to be the end of the pore. This keyword can be used to specify the radius above which the program regards a result as indicating that the end of the pore has been reached. The default value is 22.0 Å. This may need to be increased for large channels or reduced for small.
  • shorto (int, optional) – Determines the output of output in the logfile; for automated processing this must be < 3. The default is 0, which shows all output. 0: Full text output, 1: All text output given except “run in progress” (i.e., detailed contemporary description of what HOLE is doing). 2: Ditto plus no graph type output - only leaving minimum radius and conductance calculations. 3: All text output other than input card mirroring and error messages turned off.
  • ignore_residues (array_like, optional) – sequence of three-letter residues that are not taken into account during the calculation; wildcards are not supported. The default is ["SOL","WAT", "TIP", "HOH", "K ", "NA ", "CL "].
  • radius (string, optional) – path to the file specifying van der Waals radii for each atom. If set to None (the default) then a set of default radii, SIMPLE2_RAD, is used (an extension of simple.rad from the HOLE distribution)
  • executable (string, optional) – Path to the hole executable (e.g. ~/hole2/exe/hole); the other programs sph_process and sos_triangle are assumed to live in the same directory as hole. If hole is found on the PATH then the bare executable name is sufficient. The default is “hole”.
  • raseed (int, optional) – integer number to start the random number generator; by default, hole will use the time of the day (default value None) but for reproducible runs (e.g., for testing) set it to an integer.

Notes

  • An alternative way to load in multiple files is a direct read from a CHARMM binary dynamics DCD coordinate file - using the dcd keyword or use HOLEtraj.
  • HOLE is very picky and does not read all DCD-like formats [1]. If in doubt, look into the logfile for error diagnostics.
profiles

After running HOLE.collect(), this dict contains all the HOLE profiles, indexed by the frame number. If only a single frame was analyzed then this will be HOLE.profiles[0]. The entries are stored in order of calculation, but one can also sort it by the keys.

Note

Duplicate keys are not possible. The last key overwrites previous values. This is arguably a bug.

HOLE_MAX_LENGTH = 70

Maximum number of characters in a filename (limitation of HOLE)

check_and_fix_long_filename(filename, tmpdir='.')[source]

Return a file name suitable for HOLE.

HOLE is limited to filenames <= HOLE.HOLE_MAX_LENGTH. This method

  1. returns filename if HOLE can process it
  2. returns a relative path (see os.path.relpath()) if that shortens the path sufficiently
  3. creates a symlink to filename (os.symlink()) in a safe temporary directory and returns the path of the symlink. The temporary directory and the symlink are stored in HOLE.tempfiles and HOLE.tempdirs and deleted when the HOLE instance is deleted or garbage collected.
Parameters:
  • filename (string) – file name to be processed
  • tmpdir (string, optional) – By default the temporary directory is created inside the current directory in order to keep that path name short. This can be changed with the tmpdir keyword (e.g. one can use “/tmp”). The default is the current directory os.path.curdir.
Returns:

path to the file that has a length less than HOLE.HOLE_MAX_LENGTH

Return type:

string

Raises:

RuntimeError – If none of the tricks for filename shortening worked. In this case, manually rename the file or recompile your version of HOLE.

collect(**kwargs)[source]

Parse the output from a HOLE run into numpy recarrays.

It can process outputs containing multiple frames (when a DCD was supplied to hole).

Output format:

frame  rxncoord   radius

The method saves the result as HOLE.profiles, a dictionary indexed by the frame number. Each entry is a numpy.recarray.

If the keyword outdir is supplied (e.g. “.”) then each profile is saved to a gzipped data file.

Parameters:
  • run (int or string, optional) – identifier, free form; default is 1
  • outdir (string, optional) – save output data under outdir/run if set to any other value but None; the default is None.
create_vmd_surface(filename='hole.vmd', **kwargs)[source]

Process HOLE output to create a smooth pore surface suitable for VMD.

Takes the sphpdb file and feeds it to sph_process and sos_triangle as described under Visualization of HOLE results.

Load the output file filename into VMD by issuing in the tcl console

source hole.vmd

The level of detail is determined by HOLE.dotden (which can be overriden by keyword dotden).

The surface will be colored so that parts that are inaccessible to water (pore radius < 1.15 Å) are red, water accessible parts (1.15 Å > pore radius < 2.30 Å) are green and wide areas (pore radius > 2.30 Å are blue).

default_ignore_residues = ['SOL', 'WAT', 'TIP', 'HOH', 'K ', 'NA ', 'CL ']

List of residues that are ignore by default. Can be changed with the ignore_residues keyword.

min_radius()

Return the minimum radius over all profiles as a function of q

plot(**kwargs)

Plot HOLE profiles \(R(\zeta)\) in a 1D graph.

Lines are colored according to the color map cmap. One graph is plotted for each trajectory frame (unless step is changed).

Parameters:
  • step (integer, optional) – only plot every step profiles; default is 1.
  • color (string or iterable, optional) – color or iterable of colors to specify graph colors; The default None will use the color map cmap instead.
  • cmap (matplotlib.colors.Colormap) – Pick colors from the matplotlib color map cmap; the default is viridis.
  • linestyle (string or iterable, optional) – linestyle supported by matplotlib.pyplot.plot(); an iterable that is not a string (e.g., ('-', '--', '.')) will be applied in turn. The default is ‘-‘ to draw a simple line.
  • yshift (float, optional) – displace each \(R(\zeta)\) profile by yshift in the \(y\)-direction for clearer visualization. The default is 0, i.e., not to shift any graph.
  • frames (integer or array_like, optional) – only plot these specific frame(s); the default None is to plot everything (see also step)
  • label (bool or string, optional) – If it is set to “_nolegend_” or False then no legend is displayed. Any other value is ignored and the frame number is used as a label. The default is True to plot the legend even though this can become rather messy.
  • ax (matplotlib.axes.Axes) – If no ax is supplied or set to None then the plot will be added to the current active axes.
  • kwargs (**kwargs) – All other kwargs are passed to matplotlib.pyplot.plot().
Returns:

ax – Axes with the plot, either ax or the current axes.

Return type:

Axes

Notes

Changed in version 0.16.0: Returns ax.

plot3D(**kwargs)

Stacked 3D graph of profiles \(R(\zeta)\).

Lines are coloured according to the colour map cmap.

Parameters:
  • step (int, optional) – only plot every step profile; default: 1
  • cmap (matplotlib.colors.Colormap, optional) – matplotlib color map ; the default is viridis
  • rmax (float, optional) – only display radii up to rmax; the default is None, which means to show all radii.
  • ylabel (string, optional) – label of the reaction coordinate axis; the default is “frames”
  • ax (matplotlib.axes.Axes) – axes instance to plot into, which must have been generated with the projection=’3d’ keyword set (see mpl_toolkits.mplot3d.axes3d and the mplot3d tutorial for details)
Returns:

ax – Axes with the plot, either ax or the current axes.

Return type:

Axes

Notes

Based on StackOverflow post 3d plots using matplotlib. Masking of the data above rmax implements the StackOverflow hack How do I set a maximum value for the z axis.

Changed in version 0.16.0: Returns ax.

run(**kwargs)[source]

Run HOLE on the input file.

sorted_profiles_iter()

Return an iterator over profiles sorted by frame/order parameter q.

The iterator produces tuples (q, profile).

class MDAnalysis.analysis.hole.HOLEtraj(universe, **kwargs)[source]

Analyze all frames in a trajectory.

The HOLE class provides a direct interface to HOLE. HOLE itself has limited support for analysing trajectories but cannot deal with all the trajectory formats understood by MDAnalysis. This class can take any universe and feed it to HOLE. It sequentially creates a temporary PDB for each frame and runs HOLE on the frame.

The trajectory can be sliced by passing the start, stop, and step keywords to HOLEtraj.run(). (hole is not fast so slicing a trajectory is recommended.)

Frames of the trajectory can be associated with order parameters (e.g., RMSD) in order to group the HOLE profiles by order parameter (see the orderparameters keyword).

Set up the HOLE analysis over a trajectory.

Parameters:
  • universe (Universe) – The input trajectory is taken from a Universe. The trajectory is converted to a sequence of PDB files and HOLE is run on each individual file.
  • orderparameters (array_like or string, optional) – Sequence or text file containing order parameters (float numbers) corresponding to the frames in the trajectory.
  • select (string, optional) – selection string for select_atoms() to select the group of atoms that is to be analysed by HOLE. The default is “protein” to include all protein residues.
  • cpoint (bool or array_like, optional) –

    Point inside the pore to start the HOLE search procedure. The default is None to select HOLE’s internal search procedure.

    If set to True then cpoint is guessed as the center_of_geometry() of the select from the first frame of the trajectory.

    If cpoint is not set or set to None then HOLE guesses it with its own algorithm (for each individual frame).

  • kwargs (**kwargs) – All other keywords are passed on to HOLE (see there for description).

Changed in version 1.0.0: Support for the start, stop, and step keywords has been removed. These should instead be passed to HOLEtraj.run().

Changed in version 1.0.0: Changed selection keyword to select

profiles

After running HOLE.collect(), this dict contains all the HOLE profiles, indexed by the frame number or the order parameter (if orderparameters was supplied). The entries are stored in order of calculation, but one can also sort it by the keys.

Note

Duplicate keys are not possible. The last key overwrites previous values. This is arguably a bug.

guess_cpoint(select='protein', **kwargs)[source]

Guess a point inside the pore.

This method simply uses the center of geometry of the selection as a guess. select is “protein” by default.

Changed in version 1.0.0: Changed selection keyword to select

min_radius()

Return the minimum radius over all profiles as a function of q

plot(**kwargs)

Plot HOLE profiles \(R(\zeta)\) in a 1D graph.

Lines are colored according to the color map cmap. One graph is plotted for each trajectory frame (unless step is changed).

Parameters:
  • step (integer, optional) – only plot every step profiles; default is 1.
  • color (string or iterable, optional) – color or iterable of colors to specify graph colors; The default None will use the color map cmap instead.
  • cmap (matplotlib.colors.Colormap) – Pick colors from the matplotlib color map cmap; the default is viridis.
  • linestyle (string or iterable, optional) – linestyle supported by matplotlib.pyplot.plot(); an iterable that is not a string (e.g., ('-', '--', '.')) will be applied in turn. The default is ‘-‘ to draw a simple line.
  • yshift (float, optional) – displace each \(R(\zeta)\) profile by yshift in the \(y\)-direction for clearer visualization. The default is 0, i.e., not to shift any graph.
  • frames (integer or array_like, optional) – only plot these specific frame(s); the default None is to plot everything (see also step)
  • label (bool or string, optional) – If it is set to “_nolegend_” or False then no legend is displayed. Any other value is ignored and the frame number is used as a label. The default is True to plot the legend even though this can become rather messy.
  • ax (matplotlib.axes.Axes) – If no ax is supplied or set to None then the plot will be added to the current active axes.
  • kwargs (**kwargs) – All other kwargs are passed to matplotlib.pyplot.plot().
Returns:

ax – Axes with the plot, either ax or the current axes.

Return type:

Axes

Notes

Changed in version 0.16.0: Returns ax.

plot3D(**kwargs)

Stacked 3D graph of profiles \(R(\zeta)\).

Lines are coloured according to the colour map cmap.

Parameters:
  • step (int, optional) – only plot every step profile; default: 1
  • cmap (matplotlib.colors.Colormap, optional) – matplotlib color map ; the default is viridis
  • rmax (float, optional) – only display radii up to rmax; the default is None, which means to show all radii.
  • ylabel (string, optional) – label of the reaction coordinate axis; the default is “frames”
  • ax (matplotlib.axes.Axes) – axes instance to plot into, which must have been generated with the projection=’3d’ keyword set (see mpl_toolkits.mplot3d.axes3d and the mplot3d tutorial for details)
Returns:

ax – Axes with the plot, either ax or the current axes.

Return type:

Axes

Notes

Based on StackOverflow post 3d plots using matplotlib. Masking of the data above rmax implements the StackOverflow hack How do I set a maximum value for the z axis.

Changed in version 0.16.0: Returns ax.

run(start=None, stop=None, step=None)[source]

Run HOLE on the whole trajectory and collect profiles.

Parameters:
  • start (int, optional) – First frame of trajectory to analyse, Default: None (first frame)
  • stop (int, optional) – Last frame of trajectory to analyse, Default: None (last frame)
  • step (int, optional) – Step between frames to analyse, Default: None (use every frame)
  • versionchanged: (.) – 1.0.0: Undocumented support for setting HOLE parameters via HOLEtraj.run() has been removed. This should now exclusively be done on HOLEtraj construction.
run_hole(pdbfile, **kwargs)[source]

Run HOLE on a single PDB file pdbfile.

sorted_profiles_iter()

Return an iterator over profiles sorted by frame/order parameter q.

The iterator produces tuples (q, profile).

4.4.1.4. Utilities

MDAnalysis.analysis.hole.write_simplerad2(filename='simple2.rad')[source]

Write the built-in radii in SIMPLE2_RAD to filename.

Does nothing if filename already exists.

Parameters:filename (string, optional) – output file name; the default is “simple2.rad”
Returns:filename – returns the name of the data file
Return type:string
MDAnalysis.analysis.hole.SIMPLE2_RAD = '\nremark: Time-stamp: <2005-11-21 13:57:55 oliver> [OB]\nremark: van der Waals radii: AMBER united atom\nremark: from Weiner et al. (1984), JACS, vol 106 pp765-768\nremark: Simple - Only use one value for each element C O H etc.\nremark: van der Waals radii\nremark: general last\nVDWR C??? ??? 1.85\nVDWR O??? ??? 1.65\nVDWR S??? ??? 2.00\nVDWR N??? ??? 1.75\nVDWR H??? ??? 1.00\nVDWR H? ??? 1.00\nVDWR P??? ??? 2.10\nremark: ASN, GLN polar H (odd names for these atoms in xplor)\nVDWR E2? GLN 1.00\nVDWR D2? ASN 1.00\nremark: amber lone pairs on sulphurs\nVDWR LP?? ??? 0.00\nremark: for some funny reason it wants radius for K even though\nremark: it is on the exclude list\nremark: Use Pauling hydration radius (Hille 2001) [OB]\nVDWR K? ??? 1.33\nVDWR NA? ??? 0.95\nVDWR CL? ??? 1.81\nremark: funny hydrogens in gA structure [OB]\nVDWR 1H? ??? 1.00\nVDWR 2H? ??? 1.00\nVDWR 3H? ??? 1.00\nremark: need bond rad for molqpt option\nBOND C??? 0.85\nBOND N??? 0.75\nBOND O??? 0.7\nBOND S??? 1.1\nBOND H??? 0.5\nBOND P??? 1.0\nBOND ???? 0.85\n'

van der Waals radii are AMBER united atom from Weiner et al. (1984), JACS, vol 106 pp765-768. Simple - Only use one value for each element C O H etc. Added radii for K+, NA+, CL- (Pauling hydration radius from Hille 2002). The data file can be written with the convenience function write_simplerad2().

Type:Built-in HOLE radii (based on simple.rad from the HOLE distribution)
MDAnalysis.analysis.hole.seq2str(v)[source]

Return sequence as a string of numbers with spaces as separators.

In the special case of None, the empty string “” is returned.

Parameters:v (sequence or array_like) –
Returns:
Return type:string
class MDAnalysis.analysis.hole.ApplicationError[source]

Raised when an external application failed.

The error code is specific for the application.

New in version 0.7.7.

References

[Smart1993]O.S. Smart, J.M. Goodfellow and B.A. Wallace. The Pore Dimensions of Gramicidin A. Biophysical Journal 65:2455-2460, 1993. DOI: 10.1016/S0006-3495(93)81293-1
[Smart1996]O.S. Smart, J.G. Neduvelil, X. Wang, B.A. Wallace, and M.S.P. Sansom. HOLE: A program for the analysis of the pore dimensions of ion channel structural models. J.Mol.Graph., 14:354–360, 1996. DOI: 10.1016/S0263-7855(97)00009-X URL http://www.holeprogram.org/
[Stelzl2014]L. S. Stelzl, P. W. Fowler, M. S. P. Sansom, and O. Beckstein. Flexible gates generate occluded intermediates in the transport cycle of LacY. J Mol Biol, 426:735–751, 2014. DOI: 10.1016/j.jmb.2013.10.024
[1](1, 2) PDB files are not the only files that hole can read. In principle, it is also able to read CHARMM DCD trajectories and generate a hole profile for each frame. However, native support for DCD in hole is patchy and not every DCD is recognized. In particular, At the moment, DCDs generated with MDAnalysis are not accepted by HOLE. If DCDs do not work, use HOLEtraj, which converts any of the any trajectory that MDAnalysis can read into a sequence of PDB files and runs HOLE on each of them in turn.