1
2
3
4
5
6 """
7 RMSD calculation
8 ================
9
10 Calculation of the root mean square distance of a protein structure
11 over the course of a MD simulation.
12
13
14 Plugin class
15 ------------
16
17 .. autoclass:: RMSD
18 :members: worker_class
19 :undoc-members:
20
21 Worker class
22 ------------
23
24 The worker class performs the analysis.
25
26 .. autoclass:: _RMSD
27 :members:
28
29
30 """
31 from __future__ import with_statement
32
33 __docformat__ = "restructuredtext en"
34
35 import os.path
36 import warnings
37
38 import gromacs
39 from gromacs.utilities import AttributeDict
40 from gromacs.analysis.core import Worker, Plugin
41
42 import logging
43 logger = logging.getLogger('gromacs.analysis.plugins.rmsd')
44
45
46
47
48
50 """RMSD worker class."""
51
53 """Set up RMSD analysis.
54
55 This is the worker class; this is where all the real analysis is done.
56 """
57 super(_RMSD, self).__init__(**kwargs)
58 if not self.simulation is None:
59 self._register_hook()
60
62 """Run when registering; requires simulation."""
63
64 super(_RMSD, self)._register_hook(**kwargs)
65 assert not self.simulation is None
66
67 self.parameters.filenames = {
68 'RMSD': self.plugindir('rmsd.xvg'),
69 }
70
71 self.parameters.figname = self.figdir('rmsd')
72
73 - def run(self, force=False, group='C-alpha', **gmxargs):
74 """Analyze trajectory and write RMSD files.
75
76 Each frame is fit to the structure in the tpr, based on the
77 C-alpha atoms. The RMSd is alculated for the supplied index
78 group *group*.
79
80 :Arguments:
81 - *group*: index group for RMSD calculation (eg C-alpha or Protein)
82 - *force*: do analysis and overwrite existing files
83 - *gmxargs*: additional keyword arguments for :func:`gromacs.g_rms` (e.g. res=True)
84 """
85 if not self.check_file_exists(self.parameters.filenames['RMSD'], resolve='warning') or force:
86 logger.info("Analyzing RMSD...")
87 gromacs.g_rms(s=self.simulation.tpr, f=self.simulation.xtc, fit="rot+trans",
88 o=self.parameters.filenames['RMSD'],
89 input=['C-alpha', group], **gmxargs)
90
92 """Collect output xvg files as :class:`gromacs.formats.XVG` objects.
93
94 :Returns: a dictionary of the results and also sets ``self.results``.
95 """
96 from gromacs.formats import XVG
97
98 logger.info("Preparing RMSD graphs as XVG objects.")
99 results = AttributeDict(RMSD=XVG(self.parameters.filenames['RMSD']))
100 self.results = results
101 return results
102
103 - def plot(self, **kwargs):
104 """Plot all results in one graph, labelled by the result keys.
105
106 :Keywords:
107 figure
108 - ``True``: save figures in the given formats
109 - "name.ext": save figure under this filename (``ext`` -> format)
110 - ``False``: only show on screen
111 formats : sequence
112 sequence of all formats that should be saved [('png', 'pdf')]
113 plotargs
114 keyword arguments for pylab.plot()
115 """
116
117 import pylab
118 figure = kwargs.pop('figure', False)
119 extensions = kwargs.pop('formats', ('pdf','png'))
120 for name,result in self.results.items():
121 kwargs['label'] = name
122 try:
123 result.plot(**kwargs)
124 except AttributeError:
125 warnings.warn("Sorry, plotting of result %(name)r is not implemented" % vars(),
126 category=UserWarning)
127 pylab.legend(loc='best')
128 if figure is True:
129 for ext in extensions:
130 self.savefig(ext=ext)
131 elif figure:
132 self.savefig(filename=figure)
133
134
135
136
137
138
139
141 """*RMSD* plugin.
142
143 Calculation of the root mean square distance (RMSD) of a protein
144 structure over the course of a MD simulation.
145
146 The trajectory is always fitted to the reference structure in the
147 tpr file.
148
149 .. class:: RMSD([name[, simulation]])
150
151 :Arguments:
152 *name* : string
153 plugin name (used to access it)
154 *simulation* : instance
155 The :class:`gromacs.analysis.Simulation` instance that owns the plugin.
156
157 """
158 worker_class = _RMSD
159