4.4.1. HOLE analysis — MDAnalysis.analysis.hole2
¶
Author: | Lily Wang |
---|---|
Year: | 2020 |
Copyright: | GNU Public License v3 |
New in version 1.0.0.
This module provides an updated interface for the HOLE suite of tools [Smart1993]
[Smart1996] to analyse an ion channel pore or transporter pathway [Stelzl2014]
as a function of time or arbitrary order parameters. It replaces MDAnalysis.analysis.hole
.
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. Module¶
4.4.1.1.1. HOLE Analysis — MDAnalysis.analysis.hole2.hole
¶
Author: | Lily Wang |
---|---|
Year: | 2020 |
Copyright: | GNU Public License v3 |
New in version 1.0.0.
This module contains the tools to interface with HOLE [Smart1993] [Smart1996] to analyse an ion channel pore or transporter pathway [Stelzl2014] .
4.4.1.1.1.1. Using HOLE on a PDB file¶
Use the :func:hole
function to run HOLE on a single PDB file. For example,
the code below runs the HOLE program installed at ‘~/hole2/exe/hole’
from MDAnalysis.tests.datafiles import PDB_HOLE
from MDAnalysis.analysis import hole2
profiles = hole2.hole(PDB_HOLE, executable='~/hole2/exe/hole')
# to create a VMD surface of the pore
hole2.create_vmd_surface(filename='hole.vmd')
profiles
is a dictionary of HOLE profiles, indexed by the frame number. If only
a PDB file is passed to the function, there will only be one profile at frame 0.
You can visualise the pore by loading your PDB file into VMD, and in
Extensions > Tk Console, type:
source hole.vmd
You can also pass a DCD trajectory with the same atoms in the same order as
your PDB file with the dcd
keyword argument. In that case, profiles
will
contain multiple HOLE profiles, indexed by frame.
The HOLE program will create some output files:
- an output file (default name: hole.out)
- an sphpdb file (default name: hole.sph)
- a file of van der Waals’ radii (if not specified with
vdwradii_file
. Default name: simple2.rad)- a symlink of your PDB or DCD files (if the original name is too long)
- the input text (if you specify
infile
)
By default (keep_files=True), these files are kept. If you would like to delete the files after the function is wrong, set keep_files=False. Keep in mind that if you delete the sphpdb file, you cannot then create a VMD surface.
4.4.1.1.1.2. Using HOLE on a trajectory¶
You can also run HOLE on a trajectory through the HoleAnalysis
class.
This behaves similarly to the hole
function, although arguments such as cpoint
and cvect
become runtime arguments for the run()
function.
The class can be set-up and run like a normal MDAnalysis analysis class:
import MDAnalysis as mda
from MDAnalysis.tests.datafiles import MULTIPDB_HOLE
from MDAnalysis.analysis import hole2
ha = hole2.HoleAnalysis(u, executable='~/hole2/exe/hole') as h2:
ha.run()
ha.create_vmd_surface(filename='hole.vmd')
The VMD surface created by the class updates the pore for each frame of the trajectory. Use it as normal by loading your trajectory in VMD and sourcing the file in the Tk Console.
You can access the actual profiles generated in the results
attribute:
print(ha.results.profiles)
Again, HOLE writes out files for each frame. If you would like to delete these files
after the analysis, you can call delete_temporary_files()
:
ha.delete_temporary_files()
Alternatively, you can use HoleAnalysis as a context manager that deletes temporary files when you are finished with the context manager:
import MDAnalysis as mda
from MDAnalysis.tests.datafiles import MULTIPDB_HOLE
from MDAnalysis.analysis import hole2
with hole2.HoleAnalysis(u, executable='~/hole2/exe/hole') as h2:
h2.run()
h2.create_vmd_surface()
4.4.1.1.1.3. Functions and classes¶
-
MDAnalysis.analysis.hole2.hole.
hole
(pdbfile, infile_text=None, infile=None, outfile='hole.out', sphpdb_file='hole.sph', vdwradii_file=None, executable='hole', tmpdir='.', sample=0.2, end_radius=22.0, cpoint=None, cvect=None, random_seed=None, ignore_residues=['SOL', 'WAT', 'TIP', 'HOH', 'K ', 'NA ', 'CL '], output_level=0, dcd=None, dcd_iniskip=0, dcd_step=1, keep_files=True)[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.
Parameters: - pdbfile (str) – The filename is used as input for HOLE in the “COORD” card of the input file. It specifies the name of a PDB coordinate file to be used. This must be in Brookhaven protein databank format or something closely approximating this. Both ATOM and HETATM records are read.
- infile_text (str, optional) – HOLE input text or template. If set to
None
, the function will create the input text from the other parameters. - infile (str, optional) – File to write the HOLE input text for later inspection. If set to
None
, the input text is not written out. - outfile (str, optional) – file name of the file collecting HOLE’s output (which can be
parsed using
collect_hole(outfile)()
. - sphpdb_file (str, optional) – path to the HOLE sph file, a PDB-like file containing the coordinates of the pore centers. The coordinates 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. 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. .sph files can be used to produce molecular graphical output from a hole run, by using the sph_process program to read the .sph file.
- vdwradii_file (str, optional) – path to the file specifying van der Waals radii for each atom. If
set to
None
, then a set of default radii,SIMPLE2_RAD
, is used (an extension ofsimple.rad
from the HOLE distribution). - executable (str, optional) – Path to the hole executable.
(e.g.
~/hole2/exe/hole
). If hole is found on thePATH
, then the bare executable name is sufficient. - tmpdir (str, optional) – The temporary directory that files can be symlinked to, to shorten the path name. HOLE can only read filenames up to a certain length.
- sample (float, optional) – distance of sample points in Å. 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.
- end_radius (float, optional) – Radius in Å, 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. This may need to be increased for large channels, or reduced for small channels.
- cpoint (array_like, 'center_of_geometry' or None, optional) – coordinates of a point inside the pore, e.g.
[12.3, 0.7, 18.55]
. If set toNone
(the default) then HOLE’s own search algorithm is used.cpoint
specifies a point which lies within the channel. For simple channels (e.g. 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 coordinates. Or if the channel structure contains water molecules or counter ions then take the coordinates of one of these (and use theignore_residues
keyword to ignore them in the pore radius calculation). If this card is not specified, then HOLE (from version 2.2) attempts to 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 far from infallible — results should be carefully 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. If this keyword isNone
(the default), then HOLE attempts to guess where the channel will be. The procedure assumes that 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. - random_seed (int, optional) – integer number to start the random number generator.
By default,
hole will use the time of the day.
For reproducible runs (e.g., for testing) set
random_seed
to an integer. - ignore_residues (array_like, optional) – sequence of three-letter residues that are not taken into account during the calculation; wildcards are not supported. Note that all residues must have 3 letters. Pad with space on the right-hand side if necessary.
- output_level (int, optional) – Determines the output of output in the
outfile
. For automated processing, this must be < 3. 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. - dcd (str, optional) – File name of CHARMM-style DCD trajectory (must be supplied together with a
matching PDB file filename) and then HOLE runs its analysis on
each frame. HOLE can not read DCD trajectories written by MDAnalysis,
which are NAMD-style (see Notes). 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
dcd_step
keyword and/or settingdcd_iniskip
to the number of frames to be skipped initially. - dcd_step (int, optional) – step size for going through the trajectory (skips
dcd_step-1
frames). - keep_files (bool, optional) – Whether to keep the HOLE output files and possible temporary symlinks after running the function.
Returns: A dictionary of
numpy.recarray
s, indexed by frame.Return type: Notes
- HOLE is very picky and does not read all DCD-like formats [1]. If in doubt, look into the outfile for error diagnostics.
New in version 1.0.
-
class
MDAnalysis.analysis.hole2.hole.
HoleAnalysis
(universe, select='protein', verbose=False, ignore_residues=['SOL', 'WAT', 'TIP', 'HOH', 'K ', 'NA ', 'CL '], vdwradii_file=None, executable='hole', sos_triangle='sos_triangle', sph_process='sph_process', tmpdir='.', cpoint=None, cvect=None, sample=0.2, end_radius=22, output_level=0, prefix=None, write_input_files=False)[source]¶ Run hole on a 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.
This class creates temporary PDB files for each frame and runs HOLE on the frame. It can be used normally, or as a context manager. If used as a context manager, the class will try to delete any temporary files created by HOLE, e.g. sphpdb files and logfiles.
with hole2.HoleAnalysis(u, executable='~/hole2/exe/hole') as h2: h2.run() h2.create_vmd_surface()
Parameters: - universe (Universe or AtomGroup) – The Universe or AtomGroup to apply the analysis to.
- select (string, optional) – The selection string to create an atom selection that the HOLE analysis is applied to.
- vdwradii_file (str, optional) – path to the file specifying van der Waals radii for each atom. If
set to
None
, then a set of default radii,SIMPLE2_RAD
, is used (an extension ofsimple.rad
from the HOLE distribution). - executable (str, optional) – Path to the hole executable.
(e.g.
~/hole2/exe/hole
). If hole is found on thePATH
, then the bare executable name is sufficient. - tmpdir (str, optional) – The temporary directory that files can be symlinked to, to shorten the path name. HOLE can only read filenames up to a certain length.
- cpoint (array_like, 'center_of_geometry' or None, optional) – coordinates of a point inside the pore, e.g.
[12.3, 0.7, 18.55]
. If set toNone
(the default) then HOLE’s own search algorithm is used.cpoint
specifies a point which lies within the channel. For simple channels (e.g. 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 coordinates. Or if the channel structure contains water molecules or counter ions then take the coordinates of one of these (and use theignore_residues
keyword to ignore them in the pore radius calculation). If this card is not specified, then HOLE (from version 2.2) attempts to 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 far from infallible — results should be carefully 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. If this keyword isNone
(the default), then HOLE attempts to guess where the channel will be. The procedure assumes that 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 Å. 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.
- end_radius (float, optional) – Radius in Å, 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. This may need to be increased for large channels, or reduced for small channels.
- output_level (int, optional) – Determines the output of output in the
outfile
. For automated processing, this must be < 3. 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. Note that all residues must have 3 letters. Pad with space on the right-hand side if necessary.
- prefix (str, optional) – Prefix for HOLE output files.
- write_input_files (bool, optional) – Whether to write out the input HOLE text as files. Files are called hole.inp.
-
results.
sphpdbs
¶ Array of sphpdb filenames
New in version 2.0.0.
Type: numpy.ndarray
-
results.
outfiles
¶ Arrau of output filenames
New in version 2.0.0.
Type: numpy.ndarray
-
results.
profiles
¶ Profiles generated by HOLE2. A dictionary of
numpy.recarray
s, indexed by frame.New in version 2.0.0.
Type: dict
-
sphpdbs
¶ Alias of
results.sphpdbs
Deprecated since version 2.0.0: This will be removed in MDAnalysis 3.0.0. Please use
results.sphpdbs
instead.Type: numpy.ndarray
-
outfiles
¶ Alias of
results.outfiles
Deprecated since version 2.0.0: This will be removed in MDAnalysis 3.0.0. Please use
results.outfiles
instead.Type: numpy.ndarray
-
profiles
¶ Alias of
results.profiles
Deprecated since version 2.0.0: This will be removed in MDAnalysis 3.0.0. Please use
results.profiles
instead.Type: dict
-
.. versionadded:: 1.0
-
.. versionchanged:: 2.0.0
sphpdbs
,outfiles
andprofiles ` are now stored in a :class:`MDAnalysis.analysis.base.Results
instance.
-
bin_radii
(frames=None, bins=100, range=None)[source]¶ Collects the pore radii into bins by reaction coordinate.
Parameters: - frames (int or iterable of ints, optional) – Profiles to include by frame. If
None
, includes all frames. - bins (int or iterable of edges, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform bin widths.
- range ((float, float), optional) – The lower and upper range of the bins.
If not provided,
range
is simply(a.min(), a.max())
, wherea
is the array of reaction coordinates. Values outside the range are ignored. The first element of the range must be less than or equal to the second.
Returns: - list of arrays of floats – List of radii present in each bin
- array of (float, float) – Edges of each bin
- frames (int or iterable of ints, optional) – Profiles to include by frame. If
-
create_vmd_surface
(filename='hole.vmd', dot_density=15, no_water_color='red', one_water_color='green', double_water_color='blue')[source]¶ Process HOLE output to create a smooth pore surface suitable for VMD.
Takes the
sphpdb
file for each frame and feeds it to sph_process and sos_triangle as described under Visualization of HOLE results.Load the output file filename into VMD in Extensions > Tk Console
source hole.vmd
The level of detail is determined by
dot_density
. The surface will be colored byno_water_color
,one_water_color
, anddouble_water_color
. You can change these in the Tk Console:set no_water_color blue
Parameters: - filename (str, optional) – file to write the pore surfaces to.
- dot_density (int, optional) – density of facets for generating a 3D pore representation.
The number controls the density of dots that will be used.
A sphere of dots is placed on each centre determined in the
Monte Carlo procedure. The actual number of dots written is
controlled by
dot_density
and thesample
level of the original analysis.dot_density
should be set between 5 (few dots per sphere) and 35 (many dots per sphere). - no_water_color (str, optional) – Color of the surface where the pore radius is too tight for a water molecule.
- one_water_color (str, optional) – Color of the surface where the pore can fit one water molecule.
- double_water_color (str, optional) – Color of the surface where the radius is at least double the minimum radius for one water molecule.
Returns: filename
with the pore surfaces.Return type:
-
gather
(frames=None, flat=False)[source]¶ Gather the fields of each profile recarray together.
Parameters: Returns: dictionary of fields
Return type:
-
guess_cpoint
()[source]¶ Guess a point inside the pore.
This method simply uses the center of geometry of the selection as a guess.
Returns: center of geometry of selected AtomGroup Return type: float
-
histogram_radii
(aggregator=<function mean>, frames=None, bins=100, range=None)[source]¶ Histograms the pore radii into bins by reaction coordinate, aggregate the radii with an aggregator function, and returns the aggregated radii and bin edges.
Parameters: - aggregator (callable, optional) – this function must take an iterable of floats and return a single value.
- frames (int or iterable of ints, optional) – Profiles to include by frame. If
None
, includes all frames. - bins (int or iterable of edges, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform bin widths.
- range ((float, float), optional) – The lower and upper range of the bins.
If not provided,
range
is simply(a.min(), a.max())
, wherea
is the array of reaction coordinates. Values outside the range are ignored. The first element of the range must be less than or equal to the second.
Returns: - array of floats – histogrammed, aggregate value of radii
- array of (float, float) – Edges of each bin
-
over_order_parameters
(order_parameters, frames=None)[source]¶ Get HOLE profiles sorted over order parameters
order_parameters
.Parameters: - order_parameters (array-like or string) – Sequence or text file containing order parameters (float numbers) corresponding to the frames in the trajectory. Must be same length as trajectory.
- frames (array-like, optional) – Selected frames to return. If
None
, returns all of them.
Returns: sorted dictionary of {order_parameter:profile}
Return type:
-
plot
(frames=None, color=None, cmap='viridis', linestyle='-', y_shift=0.0, label=True, ax=None, legend_loc='best', **kwargs)[source]¶ Plot HOLE profiles \(R(\zeta)\) in a 1D graph.
Lines are colored according to the specified
color
or drawn from the color mapcmap
. One line is plotted for each trajectory frame.Parameters: - frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color or colors for the plot. If
None
, colors are drawn fromcmap
. - cmap (str, optional) – color map to make colors for the plot if
color
is not given. Names should be from thematplotlib.pyplot.cm
module. - linestyle (str or array-like, optional) – Line style for the plot.
- y_shift (float, optional) – displace each \(R(\zeta)\) profile by
y_shift
in the \(y\)-direction for clearer visualization. - label (bool or string, optional) – If
False
then no legend is displayed. - ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - legend_loc (str, optional) – Location of the legend.
- 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: - frames (array-like, optional) – Frames to plot. If
-
plot3D
(frames=None, color=None, cmap='viridis', linestyle='-', ax=None, r_max=None, ylabel='Frames', **kwargs)[source]¶ Stacked 3D graph of profiles \(R(\zeta)\).
Lines are colored according to the specified
color
or drawn from the color mapcmap
. One line is plotted for each trajectory frame.Parameters: - frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color or colors for the plot. If
None
, colors are drawn fromcmap
. - cmap (str, optional) – color map to make colors for the plot if
color
is not given. Names should be from thematplotlib.pyplot.cm
module. - linestyle (str or array-like, optional) – Line style for the plot.
- r_max (float, optional) – only display radii up to
r_max
. IfNone
, all radii are plotted. - ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - ylabel (str, optional) – Y-axis label.
- **kwargs – All other kwargs are passed to
matplotlib.pyplot.plot()
.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type: Axes3D
- frames (array-like, optional) – Frames to plot. If
-
plot3D_order_parameters
(order_parameters, frames=None, color=None, cmap='viridis', linestyle='-', ax=None, r_max=None, ylabel='Order parameter', **kwargs)[source]¶ Plot HOLE radii over order parameters as a 3D graph.
Lines are colored according to the specified
color
or drawn from the color mapcmap
. One line is plotted for each trajectory frame.Parameters: - order_parameters (array-like or string) – Sequence or text file containing order parameters(float numbers) corresponding to the frames in the trajectory. Must be same length as trajectory.
- frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color or colors for the plot. If
None
, colors are drawn fromcmap
. - cmap (str, optional) – color map to make colors for the plot if
color
is not given. Names should be from thematplotlib.pyplot.cm
module. - linestyle (str or array-like, optional) – Line style for the plot.
- ax (: class: matplotlib.axes.Axes) – If no ax is supplied or set to
None
then the plot will be added to the current active axes. - r_max (float, optional) – only display radii up to
r_max
. IfNone
, all radii are plotted. - ylabel (str, optional) – Y-axis label.
- **kwargs – All other kwargs are passed to: func: matplotlib.pyplot.plot.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type: : class: ~mpl_toolkits.mplot3d.Axes3D
-
plot_mean_profile
(bins=100, range=None, frames=None, color='blue', linestyle='-', ax=None, xlabel='Frame', fill_alpha=0.3, n_std=1, legend=True, legend_loc='best', **kwargs)[source]¶ Collects the pore radii into bins by reaction coordinate.
Parameters: - frames (int or iterable of ints, optional) – Profiles to include by frame. If
None
, includes all frames. - bins (int or iterable of edges, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform bin widths.
- range ((float, float), optional) – The lower and upper range of the bins.
If not provided,
range
is simply(a.min(), a.max())
, wherea
is the array of reaction coordinates. Values outside the range are ignored. The first element of the range must be less than or equal to the second. - color (str or array-like, optional) – Color for the plot.
- linestyle (str or array-like, optional) – Line style for the plot.
- ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - xlabel (str, optional) – X-axis label.
- fill_alpha (float, optional) – Opacity of filled standard deviation area
- n_std (int, optional) – Number of standard deviations from the mean to fill between.
- legend (bool, optional) – Whether to plot a legend.
- legend_loc (str, optional) – Location of legend.
- **kwargs – All other kwargs are passed to
matplotlib.pyplot.plot()
.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type: - frames (int or iterable of ints, optional) – Profiles to include by frame. If
-
plot_order_parameters
(order_parameters, aggregator=<built-in function min>, frames=None, color='blue', linestyle='-', ax=None, ylabel='Minimum HOLE pore radius $r$ ($\\AA$)', xlabel='Order parameter', **kwargs)[source]¶ Plot HOLE radii over order parameters. This function needs an
aggregator
function to reduce theradius
array to a single value, e.g.min
,max
, ornp.mean
.Parameters: - order_parameters (array-like or string) – Sequence or text file containing order parameters (float numbers) corresponding to the frames in the trajectory. Must be same length as trajectory.
- aggregator (callable, optional) – Function applied to the radius array of each profile to reduce it to one representative value.
- frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color for the plot.
- linestyle (str or array-like, optional) – Line style for the plot.
- ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - xlabel (str, optional) – X-axis label.
- ylabel (str, optional) – Y-axis label.
- **kwargs – All other kwargs are passed to
matplotlib.pyplot.plot()
.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type:
-
run
(start=None, stop=None, step=None, verbose=None, random_seed=None)[source]¶ Perform the calculation
Parameters: - start (int, optional) – start frame of analysis
- stop (int, optional) – stop frame of analysis
- step (int, optional) – number of frames to skip between each analysed frame
- verbose (bool, optional) – Turn on verbosity
- random_seed (int, optional) – integer number to start the random number generator.
By default,
hole will use the time of the day.
For reproducible runs (e.g., for testing) set
random_seed
to an integer.
References
[Smart1993] | (1, 2) 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] | (1, 2) 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] | (1, 2) 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. To overcome this
PDB / DCD limitation, use HoleAnalysis which creates
temporary PDB files for each frame of a
Universe or
AtomGroup and runs
:func:hole on each of them. |
-
class
MDAnalysis.analysis.hole2.hole.
HoleAnalysis
(universe, select='protein', verbose=False, ignore_residues=['SOL', 'WAT', 'TIP', 'HOH', 'K ', 'NA ', 'CL '], vdwradii_file=None, executable='hole', sos_triangle='sos_triangle', sph_process='sph_process', tmpdir='.', cpoint=None, cvect=None, sample=0.2, end_radius=22, output_level=0, prefix=None, write_input_files=False)[source] Run hole on a 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.
This class creates temporary PDB files for each frame and runs HOLE on the frame. It can be used normally, or as a context manager. If used as a context manager, the class will try to delete any temporary files created by HOLE, e.g. sphpdb files and logfiles.
with hole2.HoleAnalysis(u, executable='~/hole2/exe/hole') as h2: h2.run() h2.create_vmd_surface()
Parameters: - universe (Universe or AtomGroup) – The Universe or AtomGroup to apply the analysis to.
- select (string, optional) – The selection string to create an atom selection that the HOLE analysis is applied to.
- vdwradii_file (str, optional) – path to the file specifying van der Waals radii for each atom. If
set to
None
, then a set of default radii,SIMPLE2_RAD
, is used (an extension ofsimple.rad
from the HOLE distribution). - executable (str, optional) – Path to the hole executable.
(e.g.
~/hole2/exe/hole
). If hole is found on thePATH
, then the bare executable name is sufficient. - tmpdir (str, optional) – The temporary directory that files can be symlinked to, to shorten the path name. HOLE can only read filenames up to a certain length.
- cpoint (array_like, 'center_of_geometry' or None, optional) – coordinates of a point inside the pore, e.g.
[12.3, 0.7, 18.55]
. If set toNone
(the default) then HOLE’s own search algorithm is used.cpoint
specifies a point which lies within the channel. For simple channels (e.g. 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 coordinates. Or if the channel structure contains water molecules or counter ions then take the coordinates of one of these (and use theignore_residues
keyword to ignore them in the pore radius calculation). If this card is not specified, then HOLE (from version 2.2) attempts to 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 far from infallible — results should be carefully 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. If this keyword isNone
(the default), then HOLE attempts to guess where the channel will be. The procedure assumes that 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 Å. 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.
- end_radius (float, optional) – Radius in Å, 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. This may need to be increased for large channels, or reduced for small channels.
- output_level (int, optional) – Determines the output of output in the
outfile
. For automated processing, this must be < 3. 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. Note that all residues must have 3 letters. Pad with space on the right-hand side if necessary.
- prefix (str, optional) – Prefix for HOLE output files.
- write_input_files (bool, optional) – Whether to write out the input HOLE text as files. Files are called hole.inp.
-
results.
sphpdbs
Array of sphpdb filenames
New in version 2.0.0.
Type: numpy.ndarray
-
results.
outfiles
Arrau of output filenames
New in version 2.0.0.
Type: numpy.ndarray
-
results.
profiles
Profiles generated by HOLE2. A dictionary of
numpy.recarray
s, indexed by frame.New in version 2.0.0.
Type: dict
-
sphpdbs
Alias of
results.sphpdbs
Deprecated since version 2.0.0: This will be removed in MDAnalysis 3.0.0. Please use
results.sphpdbs
instead.Type: numpy.ndarray
-
outfiles
Alias of
results.outfiles
Deprecated since version 2.0.0: This will be removed in MDAnalysis 3.0.0. Please use
results.outfiles
instead.Type: numpy.ndarray
-
profiles
Alias of
results.profiles
Deprecated since version 2.0.0: This will be removed in MDAnalysis 3.0.0. Please use
results.profiles
instead.Type: dict
-
.. versionadded:: 1.0
-
.. versionchanged:: 2.0.0
sphpdbs
,outfiles
andprofiles ` are now stored in a :class:`MDAnalysis.analysis.base.Results
instance.
-
bin_radii
(frames=None, bins=100, range=None)[source] Collects the pore radii into bins by reaction coordinate.
Parameters: - frames (int or iterable of ints, optional) – Profiles to include by frame. If
None
, includes all frames. - bins (int or iterable of edges, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform bin widths.
- range ((float, float), optional) – The lower and upper range of the bins.
If not provided,
range
is simply(a.min(), a.max())
, wherea
is the array of reaction coordinates. Values outside the range are ignored. The first element of the range must be less than or equal to the second.
Returns: - list of arrays of floats – List of radii present in each bin
- array of (float, float) – Edges of each bin
- frames (int or iterable of ints, optional) – Profiles to include by frame. If
-
create_vmd_surface
(filename='hole.vmd', dot_density=15, no_water_color='red', one_water_color='green', double_water_color='blue')[source] Process HOLE output to create a smooth pore surface suitable for VMD.
Takes the
sphpdb
file for each frame and feeds it to sph_process and sos_triangle as described under Visualization of HOLE results.Load the output file filename into VMD in Extensions > Tk Console
source hole.vmd
The level of detail is determined by
dot_density
. The surface will be colored byno_water_color
,one_water_color
, anddouble_water_color
. You can change these in the Tk Console:set no_water_color blue
Parameters: - filename (str, optional) – file to write the pore surfaces to.
- dot_density (int, optional) – density of facets for generating a 3D pore representation.
The number controls the density of dots that will be used.
A sphere of dots is placed on each centre determined in the
Monte Carlo procedure. The actual number of dots written is
controlled by
dot_density
and thesample
level of the original analysis.dot_density
should be set between 5 (few dots per sphere) and 35 (many dots per sphere). - no_water_color (str, optional) – Color of the surface where the pore radius is too tight for a water molecule.
- one_water_color (str, optional) – Color of the surface where the pore can fit one water molecule.
- double_water_color (str, optional) – Color of the surface where the radius is at least double the minimum radius for one water molecule.
Returns: filename
with the pore surfaces.Return type:
-
delete_temporary_files
()[source] Delete temporary files
-
gather
(frames=None, flat=False)[source] Gather the fields of each profile recarray together.
Parameters: Returns: dictionary of fields
Return type:
-
guess_cpoint
()[source] Guess a point inside the pore.
This method simply uses the center of geometry of the selection as a guess.
Returns: center of geometry of selected AtomGroup Return type: float
-
histogram_radii
(aggregator=<function mean>, frames=None, bins=100, range=None)[source] Histograms the pore radii into bins by reaction coordinate, aggregate the radii with an aggregator function, and returns the aggregated radii and bin edges.
Parameters: - aggregator (callable, optional) – this function must take an iterable of floats and return a single value.
- frames (int or iterable of ints, optional) – Profiles to include by frame. If
None
, includes all frames. - bins (int or iterable of edges, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform bin widths.
- range ((float, float), optional) – The lower and upper range of the bins.
If not provided,
range
is simply(a.min(), a.max())
, wherea
is the array of reaction coordinates. Values outside the range are ignored. The first element of the range must be less than or equal to the second.
Returns: - array of floats – histogrammed, aggregate value of radii
- array of (float, float) – Edges of each bin
-
min_radius
()[source] Return the minimum radius over all profiles as a function of q
-
over_order_parameters
(order_parameters, frames=None)[source] Get HOLE profiles sorted over order parameters
order_parameters
.Parameters: - order_parameters (array-like or string) – Sequence or text file containing order parameters (float numbers) corresponding to the frames in the trajectory. Must be same length as trajectory.
- frames (array-like, optional) – Selected frames to return. If
None
, returns all of them.
Returns: sorted dictionary of {order_parameter:profile}
Return type:
-
plot
(frames=None, color=None, cmap='viridis', linestyle='-', y_shift=0.0, label=True, ax=None, legend_loc='best', **kwargs)[source] Plot HOLE profiles \(R(\zeta)\) in a 1D graph.
Lines are colored according to the specified
color
or drawn from the color mapcmap
. One line is plotted for each trajectory frame.Parameters: - frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color or colors for the plot. If
None
, colors are drawn fromcmap
. - cmap (str, optional) – color map to make colors for the plot if
color
is not given. Names should be from thematplotlib.pyplot.cm
module. - linestyle (str or array-like, optional) – Line style for the plot.
- y_shift (float, optional) – displace each \(R(\zeta)\) profile by
y_shift
in the \(y\)-direction for clearer visualization. - label (bool or string, optional) – If
False
then no legend is displayed. - ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - legend_loc (str, optional) – Location of the legend.
- 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: - frames (array-like, optional) – Frames to plot. If
-
plot3D
(frames=None, color=None, cmap='viridis', linestyle='-', ax=None, r_max=None, ylabel='Frames', **kwargs)[source] Stacked 3D graph of profiles \(R(\zeta)\).
Lines are colored according to the specified
color
or drawn from the color mapcmap
. One line is plotted for each trajectory frame.Parameters: - frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color or colors for the plot. If
None
, colors are drawn fromcmap
. - cmap (str, optional) – color map to make colors for the plot if
color
is not given. Names should be from thematplotlib.pyplot.cm
module. - linestyle (str or array-like, optional) – Line style for the plot.
- r_max (float, optional) – only display radii up to
r_max
. IfNone
, all radii are plotted. - ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - ylabel (str, optional) – Y-axis label.
- **kwargs – All other kwargs are passed to
matplotlib.pyplot.plot()
.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type: Axes3D
- frames (array-like, optional) – Frames to plot. If
-
plot3D_order_parameters
(order_parameters, frames=None, color=None, cmap='viridis', linestyle='-', ax=None, r_max=None, ylabel='Order parameter', **kwargs)[source] Plot HOLE radii over order parameters as a 3D graph.
Lines are colored according to the specified
color
or drawn from the color mapcmap
. One line is plotted for each trajectory frame.Parameters: - order_parameters (array-like or string) – Sequence or text file containing order parameters(float numbers) corresponding to the frames in the trajectory. Must be same length as trajectory.
- frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color or colors for the plot. If
None
, colors are drawn fromcmap
. - cmap (str, optional) – color map to make colors for the plot if
color
is not given. Names should be from thematplotlib.pyplot.cm
module. - linestyle (str or array-like, optional) – Line style for the plot.
- ax (: class: matplotlib.axes.Axes) – If no ax is supplied or set to
None
then the plot will be added to the current active axes. - r_max (float, optional) – only display radii up to
r_max
. IfNone
, all radii are plotted. - ylabel (str, optional) – Y-axis label.
- **kwargs – All other kwargs are passed to: func: matplotlib.pyplot.plot.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type: : class: ~mpl_toolkits.mplot3d.Axes3D
-
plot_mean_profile
(bins=100, range=None, frames=None, color='blue', linestyle='-', ax=None, xlabel='Frame', fill_alpha=0.3, n_std=1, legend=True, legend_loc='best', **kwargs)[source] Collects the pore radii into bins by reaction coordinate.
Parameters: - frames (int or iterable of ints, optional) – Profiles to include by frame. If
None
, includes all frames. - bins (int or iterable of edges, optional) – If bins is an int, it defines the number of equal-width bins in the given range. If bins is a sequence, it defines a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform bin widths.
- range ((float, float), optional) – The lower and upper range of the bins.
If not provided,
range
is simply(a.min(), a.max())
, wherea
is the array of reaction coordinates. Values outside the range are ignored. The first element of the range must be less than or equal to the second. - color (str or array-like, optional) – Color for the plot.
- linestyle (str or array-like, optional) – Line style for the plot.
- ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - xlabel (str, optional) – X-axis label.
- fill_alpha (float, optional) – Opacity of filled standard deviation area
- n_std (int, optional) – Number of standard deviations from the mean to fill between.
- legend (bool, optional) – Whether to plot a legend.
- legend_loc (str, optional) – Location of legend.
- **kwargs – All other kwargs are passed to
matplotlib.pyplot.plot()
.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type: - frames (int or iterable of ints, optional) – Profiles to include by frame. If
-
plot_order_parameters
(order_parameters, aggregator=<built-in function min>, frames=None, color='blue', linestyle='-', ax=None, ylabel='Minimum HOLE pore radius $r$ ($\\AA$)', xlabel='Order parameter', **kwargs)[source] Plot HOLE radii over order parameters. This function needs an
aggregator
function to reduce theradius
array to a single value, e.g.min
,max
, ornp.mean
.Parameters: - order_parameters (array-like or string) – Sequence or text file containing order parameters (float numbers) corresponding to the frames in the trajectory. Must be same length as trajectory.
- aggregator (callable, optional) – Function applied to the radius array of each profile to reduce it to one representative value.
- frames (array-like, optional) – Frames to plot. If
None
, plots all of them. - color (str or array-like, optional) – Color for the plot.
- linestyle (str or array-like, optional) – Line style for the plot.
- ax (
matplotlib.axes.Axes
) – If no ax is supplied or set toNone
then the plot will be added to the current active axes. - xlabel (str, optional) – X-axis label.
- ylabel (str, optional) – Y-axis label.
- **kwargs – All other kwargs are passed to
matplotlib.pyplot.plot()
.
Returns: ax – Axes with the plot, either ax or the current axes.
Return type:
-
run
(start=None, stop=None, step=None, verbose=None, random_seed=None)[source] Perform the calculation
Parameters: - start (int, optional) – start frame of analysis
- stop (int, optional) – stop frame of analysis
- step (int, optional) – number of frames to skip between each analysed frame
- verbose (bool, optional) – Turn on verbosity
- random_seed (int, optional) – integer number to start the random number generator.
By default,
hole will use the time of the day.
For reproducible runs (e.g., for testing) set
random_seed
to an integer.
-
MDAnalysis.analysis.hole2.hole.
hole
(pdbfile, infile_text=None, infile=None, outfile='hole.out', sphpdb_file='hole.sph', vdwradii_file=None, executable='hole', tmpdir='.', sample=0.2, end_radius=22.0, cpoint=None, cvect=None, random_seed=None, ignore_residues=['SOL', 'WAT', 'TIP', 'HOH', 'K ', 'NA ', 'CL '], output_level=0, dcd=None, dcd_iniskip=0, dcd_step=1, keep_files=True)[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.
Parameters: - pdbfile (str) – The filename is used as input for HOLE in the “COORD” card of the input file. It specifies the name of a PDB coordinate file to be used. This must be in Brookhaven protein databank format or something closely approximating this. Both ATOM and HETATM records are read.
- infile_text (str, optional) – HOLE input text or template. If set to
None
, the function will create the input text from the other parameters. - infile (str, optional) – File to write the HOLE input text for later inspection. If set to
None
, the input text is not written out. - outfile (str, optional) – file name of the file collecting HOLE’s output (which can be
parsed using
collect_hole(outfile)()
. - sphpdb_file (str, optional) – path to the HOLE sph file, a PDB-like file containing the coordinates of the pore centers. The coordinates 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. 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. .sph files can be used to produce molecular graphical output from a hole run, by using the sph_process program to read the .sph file.
- vdwradii_file (str, optional) – path to the file specifying van der Waals radii for each atom. If
set to
None
, then a set of default radii,SIMPLE2_RAD
, is used (an extension ofsimple.rad
from the HOLE distribution). - executable (str, optional) – Path to the hole executable.
(e.g.
~/hole2/exe/hole
). If hole is found on thePATH
, then the bare executable name is sufficient. - tmpdir (str, optional) – The temporary directory that files can be symlinked to, to shorten the path name. HOLE can only read filenames up to a certain length.
- sample (float, optional) – distance of sample points in Å. 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.
- end_radius (float, optional) – Radius in Å, 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. This may need to be increased for large channels, or reduced for small channels.
- cpoint (array_like, 'center_of_geometry' or None, optional) – coordinates of a point inside the pore, e.g.
[12.3, 0.7, 18.55]
. If set toNone
(the default) then HOLE’s own search algorithm is used.cpoint
specifies a point which lies within the channel. For simple channels (e.g. 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 coordinates. Or if the channel structure contains water molecules or counter ions then take the coordinates of one of these (and use theignore_residues
keyword to ignore them in the pore radius calculation). If this card is not specified, then HOLE (from version 2.2) attempts to 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 far from infallible — results should be carefully 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. If this keyword isNone
(the default), then HOLE attempts to guess where the channel will be. The procedure assumes that 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. - random_seed (int, optional) – integer number to start the random number generator.
By default,
hole will use the time of the day.
For reproducible runs (e.g., for testing) set
random_seed
to an integer. - ignore_residues (array_like, optional) – sequence of three-letter residues that are not taken into account during the calculation; wildcards are not supported. Note that all residues must have 3 letters. Pad with space on the right-hand side if necessary.
- output_level (int, optional) – Determines the output of output in the
outfile
. For automated processing, this must be < 3. 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. - dcd (str, optional) – File name of CHARMM-style DCD trajectory (must be supplied together with a
matching PDB file filename) and then HOLE runs its analysis on
each frame. HOLE can not read DCD trajectories written by MDAnalysis,
which are NAMD-style (see Notes). 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
dcd_step
keyword and/or settingdcd_iniskip
to the number of frames to be skipped initially. - dcd_step (int, optional) – step size for going through the trajectory (skips
dcd_step-1
frames). - keep_files (bool, optional) – Whether to keep the HOLE output files and possible temporary symlinks after running the function.
Returns: A dictionary of
numpy.recarray
s, indexed by frame.Return type: Notes
- HOLE is very picky and does not read all DCD-like formats [1]. If in doubt, look into the outfile for error diagnostics.
New in version 1.0.
4.4.1.2. Utility functions and templates¶
4.4.1.2.1. HOLE Analysis — MDAnalysis.analysis.hole2.helper
¶
Author: | Lily Wang |
---|---|
Year: | 2020 |
Copyright: | GNU Public License v3 |
New in version 1.0.
Helper functions used in MDAnalysis.analysis.hole2.hole
-
MDAnalysis.analysis.hole2.utils.
check_and_fix_long_filename
(filename, tmpdir='.', max_length=70, make_symlink=True)[source]¶ Return a file name suitable for HOLE.
HOLE is limited to filenames <=
max_length
. This method- returns filename if HOLE can process it
- returns a relative path (see
os.path.relpath()
) if that shortens the - path sufficiently
- returns a relative path (see
- creates a symlink to filename (
os.symlink()
) in a safe temporary - directory and returns the path of the symlink.
- creates a symlink to filename (
Parameters: - filename (str) – file name to be processed
- tmpdir (str, 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
max_length
Return type: Raises: RuntimeError
– If none of the tricks for filename shortening worked. In this case, manually rename the file or recompile your version of HOLE.
-
MDAnalysis.analysis.hole2.utils.
collect_hole
(outfile='hole.out')[source]¶ Collect data from HOLE output
Parameters: outfile (str, optional) – HOLE output file to read. Default: ‘hole.out’ Returns: Dictionary of HOLE profiles as record arrays Return type: dict
-
MDAnalysis.analysis.hole2.utils.
create_vmd_surface
(sphpdb='hole.sph', filename=None, sph_process='sph_process', sos_triangle='sos_triangle', dot_density=15)[source]¶ Create VMD surface file from sphpdb file.
Parameters: - sphpdb (str, optional) – sphpdb file to read. Default: ‘hole.sph’
- filename (str, optional) – output VMD surface file. If
None
, a temporary file is generated. Default:None
- sph_process (str, optional) – Executable for
sph_process
program. Default: ‘sph_process’ - sos_triangle (str, optional) – Executable for
sos_triangle
program. Default: ‘sos_triangle’ - dot_density (int, optional) – density of facets for generating a 3D pore representation.
The number controls the density of dots that will be used.
A sphere of dots is placed on each centre determined in the
Monte Carlo procedure. The actual number of dots written is
controlled by
dot_density
and thesample
level of the original analysis.dot_density
should be set between 5 (few dots per sphere) and 35 (many dots per sphere). Default: 15
Returns: the output filename for the VMD surface
Return type:
-
MDAnalysis.analysis.hole2.utils.
run_hole
(outfile, infile_text, executable)[source]¶ Run the HOLE program.
Parameters: Returns: Output file name
Return type:
-
MDAnalysis.analysis.hole2.utils.
write_simplerad2
(filename='simple2.rad')[source]¶ Write the built-in radii in
SIMPLE2_RAD
to filename.Does nothing if filename already exists.
Parameters: filename (str, optional) – output file name; the default is “simple2.rad” Returns: filename – returns the name of the data file Return type: str
4.4.1.2.2. HOLE Analysis — MDAnalysis.analysis.hole2.templates
¶
Author: | Lily Wang |
---|---|
Year: | 2020 |
Copyright: | GNU Public License v3 |
New in version 1.0.
Templates used in MDAnalysis.analysis.hole2.hole
-
MDAnalysis.analysis.hole2.templates.
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)