Package gromacs :: Package analysis :: Package plugins :: Module mindistances
[hide private]
[frames] | no frames]

Source Code for Module gromacs.analysis.plugins.mindistances

  1  # $Id$ 
  2  # Copyright (c) 2009 Oliver Beckstein <orbeckst@gmail.com> 
  3  # Released under the GNU Public License 3 (or higher, your choice) 
  4  # See the file COPYING for details. 
  5   
  6  """ 
  7  MinDistance plugin 
  8  =============== 
  9   
 10  Time series of a selected set of shortest distances. 
 11   
 12  Plugin class 
 13  ------------ 
 14   
 15  .. autoclass:: MinDistances 
 16     :members: worker_class 
 17     :undoc-members: 
 18   
 19  Worker class 
 20  ------------ 
 21   
 22  The worker class performs the analysis. 
 23   
 24  .. autoclass:: _MinDistances 
 25     :members: 
 26   
 27   
 28  """ 
 29  __docformat__ = "restructuredtext en" 
 30   
 31  import gromacs 
 32  from gromacs.utilities import asiterable 
 33  from gromacs.analysis.core import Plugin 
 34  from distances import _Distances 
 35   
 36  # Worker classes that are registered via Plugins (see below) 
 37  # ---------------------------------------------------------- 
 38  # These must be defined before the plugins. 
 39   
40 -class _MinDistances(_Distances):
41 """Analysis of (multiple) minimum distances.""" 42 43 #: list of results (not used at the moment, see _register_hook()) 44 names = ["distance", "contacts"] 45 46 default_plot_columns = Ellipsis # plot everything by default 47
48 - def _register_hook(self, **kwargs):
49 """Run when registering; requires simulation. 50 51 Defines output files (note that we overwrite the 52 parameters.filenames and figname that super might have set). 53 """ 54 55 super(_MinDistances, self)._register_hook(**kwargs) 56 assert not self.simulation is None 57 58 # output filenames for g_mindist 59 self.parameters.filenames = {'contacts': self.plugindir('contacts.xvg'), 60 'distance': self.plugindir('distance.xvg'), 61 } 62 # default filename for the combined plot 63 self.parameters.figname = self.figdir('mindist')
64
65 - def run(self,**kwargs):
66 """Run ``g_mindist `` to compute distances between multiple groups. 67 68 Additional arguments can be provided (e.g. ``-b`` or ``-e``) 69 but an error will result if one tries to set parameters that 70 are already being set by the method itself such as ``-s`` or 71 ``-d``; one must to provide the appropriate values to the 72 class constructor. 73 74 If the primary output file already exists then no data are generated 75 and the method returns immediately unless one sets *force* = ``True``. 76 """ 77 force = kwargs.pop('force',False) 78 if not force and \ 79 self.check_file_exists(self.parameters.filenames['distance'],resolve='warn'): 80 return 81 indexgroups = self.parameters.indexgroups 82 ngroups = len(indexgroups) - 1 # number of secondary groups 83 kwargs.setdefault('o', None) # set to True if default output is required, or 84 kwargs.setdefault('_or', None) # add filenames to self.parameters.filenames 85 gromacs.g_mindist(s=self.simulation.tpr, n=self.parameters.ndx, 86 f=self.simulation.xtc, d=self.parameters.cutoff, 87 od=self.parameters.filenames['distance'], 88 on=self.parameters.filenames['contacts'], 89 ng=ngroups, input=indexgroups, 90 **kwargs)
91 92 93 94 95 96 # Public classes that register the worker classes 97 #------------------------------------------------ 98
99 -class MinDistances(Plugin):
100 """*MinDistances* plugin. 101 102 The minimum distances between the members of at least two index 103 groups and the number of contacts are calculated for each time 104 step and written to files. 105 106 .. class:: Distances(groups, ndx, [cutoff, [, name[, simulation]]]) 107 108 :Arguments: 109 name : string 110 plugin name (used to access it) 111 simulation : instance 112 The :class:`gromacs.analysis.Simulation` instance that owns the plugin. 113 groups : list of index group names 114 The first entry is the *primary group*. All other entries 115 are *secondary groups* and the plugin calculates the minimum distance 116 between members of the primary group and the members of each 117 secondary group. 118 ndx : index filename or list 119 All index files that contain the listed groups. 120 cutoff : float 121 A contact is recorded if the distance is <cutoff [0.6 nm] 122 123 Example: 124 125 Generate index files with the groups of interest, for instance 126 with :class:`gromacs.cbook.IndexBuilder`:: 127 128 from gromacs.cbook import IndexBuilder 129 A_grp, A_ndx = IndexBuilder(tpr, ['@a 62549 & r NA'], names=['Na1_ion'], offset=-9, 130 out_ndx='Na1.ndx', name_all="Na1").combine() 131 B = IndexBuilder(tpr, ['S312:OG','T313:OG1','A38:O','I41:O','A309:O'], offset=-9, 132 out_ndx='Na1_site.ndx', name_all="Na1_site") 133 B_grp, B_ndx = B.combine() 134 all_ndx_files = [A_ndx, B_ndx] 135 136 To calculate the distance between "Na1" and the "Na1_site", create an instance with 137 the appropriate parameters and add them to a :class:`gromacs.analysis.Simulation` instance:: 138 139 dist_Na1_site = Distances(name='Dsite', groups=['Na1', 'Na1_site'], ndx=all_ndx_files) 140 S.add_plugin(dist_Na1_site) 141 142 To calculate the individual distances:: 143 144 dist_Na1_res = Distances(name='Dres', groups=['Na1']+B.names, ndx=all_ndx_files) 145 S.add_plugin(dist_Na1_res) 146 147 (Keeping the second IndexBuilder instance ``B`` allows us to directly 148 use all groups without typing them, ``B.names = ['A309_O', 'S312_OG', 'I41_O', 149 'T313_OG1', 'A38_O']``.) 150 151 152 """ 153 worker_class = _MinDistances
154