4.11.1.1.1. Generation and Analysis of X3DNA helicoidal parameter profiles — MDAnalysis.analysis.legacy.x3dna
- Author:
Elizabeth Denning
- Year:
2013-2014
- Copyright:
Lesser GNU Public License v2.1+
Added in version 0.8.
Changed in version 0.16.0: This module is difficult to test due to restrictions on the X3DNA code. It
is therefore considered unmaintained and legacy code. It was moved to the
MDAnalysis.analysis.legacy
package (see issue 906)
With the help of this module, X3DNA can be run on frames in a trajectory. Data can be combined and analyzed. X3DNA [1][2] must be installed separately.
References
4.11.1.1.1.1. Example applications
4.11.1.1.1.1.1. Single structure
B-DNA structure:
from MDAnalysis.analysis.x3dna import X3DNA, X3DNAtraj
from MDAnalysis.tests.datafiles import PDB_X3DNA
# set path to your x3dna binary in bashrc file
H = X3DNA(PDB_X3DNA, executable="x3dna_ensemble analyze -b 355d.bps -p pdbfile")
H.run()
H.collect()
H.plot()
4.11.1.1.1.1.2. Trajectory
Analyzing a trajectory:
u = MDAnalysis.Universe(psf, trajectory)
H = X3DNAtraj(u, ...)
H.run()
H.plot()
H.save()
The profiles are available as the attribute X3DNAtraj.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.11.1.1.1.2. Analysis classes
- class MDAnalysis.analysis.legacy.x3dna.X3DNA(**kwds)[source]
Run X3DNA on a single frame or a DCD trajectory.
Only a subset of all X3DNA control parameters is supported and can be set with keyword arguments. For further details on X3DNA see the X3DNA docs.
Running X3DNA with the
X3DNA
class is a 3-step process:set up the class with all desired parameters
run X3DNA with
X3DNA.run()
collect the data from the output file with
X3DNA.collect()
The class also provides some simple plotting functions of the collected data such as
X3DNA.plot()
orX3DNA.plot3D()
.When methods return helicoidal basepair parameter as lists, then the order is always
index
parameter
0
shear
1
stretch
2
stagger
3
buckle
4
propeller
5
opening
6
shift
7
slide
8
rise
9
tilt
10
roll
11
twist
Added in version 0.8.
Deprecated since version 2.7.0: X3DNA will be removed in 3.0.0.
__init__ is deprecated!
Set up parameters to run X3DNA on PDB filename.
- Parameters:
filename (str) – The filename is used as input for X3DNA in the xdna_ensemble command. 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.
executable (str (optional)) – Path to the xdna_ensemble executable directories (e.g.
/opt/x3dna/2.1 and /opt/x3dna/2.1/bin
) must be set and then added to export in bashrc file. See X3DNA documentation for set-up instructions.x3dna_param (bool (optional)) – Determines whether base step or base pair parameters will be calculated. If
True
(default) then stacked base step parameters will be analyzed. IfFalse
then stacked base pair parameters will be analyzed.logfile (str (optional)) – Write output from X3DNA to logfile (default: “bp_step.par”)
See also
Deprecated since version 2.7.0: X3DNA module is deprecated and will be removed inMDAnalysis 3.0.0, see #3788 __init__ will be removed in release 3.0.0.
- profiles
x3dna_ensemble analyze -b 355d.bps -p pdbfile attribute
: After runningX3DNA.collect()
, this dict contains all the X3DNA profiles, indexed by the frame number. If only a single frame was analyzed then this will beX3DNA.profiles[0]
. Note that the order is random; one needs to sort the keys first.
- collect(**kwargs)[source]
Parse the output from a X3DNA run into numpy recarrays.
Can deal with outputs containing multiple frames.
The method saves the result as
X3DNA.profiles
, a dictionary indexed by the frame number. Each entry is anp.recarray
.If the keyword outdir is supplied (e.g. “.”) then each profile is saved to a gzipped data file.
- mean()
Returns the mean value for the base parameters.
- Returns:
The list contains the means for the helicoidal parameters. The order is
[shear, stretch, stagger, buckle, propeller, opening, shift, slide, rise, tilt, roll, twist]
.- Return type:
- mean_std()
Returns the mean and standard deviation of base parameters.
- plot(**kwargs)
Plot time-averaged base parameters for each basse pair in a 1D graph.
One plot is produced for each parameter. It shows the the mean and standard deviation for each individual base pair. Each plot is saved to PNG file with name “<parameter_name>.png”.
- Parameters:
ax (matplotlib.pyplot.Axes (optional)) – Provide ax to have all plots plotted in the same axes.
- save(filename='x3dna.pickle')
Save
profiles
as a Python pickle file filename.Load profiles dictionary with
import cPickle profiles = cPickle.load(open(filename))
- sorted_profiles_iter()
Return an iterator over profiles sorted by frame/order parameter.
The iterator produces tuples
(q, profile)
. Typically, q is the frame number.
- class MDAnalysis.analysis.legacy.x3dna.X3DNAtraj(**kwds)[source]
Analyze all frames in a trajectory.
The
X3DNA
class provides a direct interface to X3DNA. X3DNA 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 X3DNA. By default it sequentially creates a PDB for each frame and runs X3DNA on the frame.Deprecated since version 2.7.0: X3DNA will be removed in 3.0.0.
__init__ is deprecated!
Set up the class.
- Parameters:
universe (Universe) – The input trajectory as part of a
Universe
; the trajectory is converted to a sequence of PDB files and X3DNA is run on each individual file. (Use the start, stop, and step keywords to slice the trajectory.)selection (str (optional)) – MDAnalysis selection string (default: “nucleic”) to select the atoms that should be analyzed.
start (int (optional))
stop (int (optional))
step (int (optional)) – frame indices to slice the trajectory as
universe.trajectory[start, stop, step]
; by default, the whole trajectory is analyzed.x3dna_param (bool (optional)) – indicates whether stacked bases or stacked base-pairs will be analyzed.
True
is bases andFalse
is stacked base-pairs [Default isTrue
].kwargs (keyword arguments (optional)) – All other keywords are passed on to
X3DNA
(see there for description).
See also
Deprecated since version 2.7.0: X3DNA module is deprecated and will be removed inMDAnalysis 3.0.0, see #3788 __init__ will be removed in release 3.0.0.
- profiles
After running
X3DNA.collect()
, this dict contains all the X3DNA profiles, indexed by the frame number.
- mean()
Returns the mean value for the base parameters.
- Returns:
The list contains the means for the helicoidal parameters. The order is
[shear, stretch, stagger, buckle, propeller, opening, shift, slide, rise, tilt, roll, twist]
.- Return type:
- mean_std()
Returns the mean and standard deviation of base parameters.
- plot(**kwargs)
Plot time-averaged base parameters for each basse pair in a 1D graph.
One plot is produced for each parameter. It shows the the mean and standard deviation for each individual base pair. Each plot is saved to PNG file with name “<parameter_name>.png”.
- Parameters:
ax (matplotlib.pyplot.Axes (optional)) – Provide ax to have all plots plotted in the same axes.
- run(**kwargs)[source]
Run X3DNA on the whole trajectory and collect profiles.
Keyword arguments start, stop, and step can be used to only analyse part of the trajectory. The defaults are the values provided to the class constructor.
- save(filename='x3dna.pickle')
Save
profiles
as a Python pickle file filename.Load profiles dictionary with
import cPickle profiles = cPickle.load(open(filename))
- sorted_profiles_iter()
Return an iterator over profiles sorted by frame/order parameter.
The iterator produces tuples
(q, profile)
. Typically, q is the frame number.