Skip to content

Commit

Permalink
removed quiet in favour of verbose
Browse files Browse the repository at this point in the history
updated usage of Analysis API, start/stop/step
now deprecated in __init__ and should be given in run method

made PCA setup respect start frame more,
if start is given then initial weights are from that frame
not frame 0
  • Loading branch information
richardjgowers committed Sep 21, 2018
1 parent 55e7475 commit 00cf3b7
Show file tree
Hide file tree
Showing 22 changed files with 239 additions and 373 deletions.
5 changes: 1 addition & 4 deletions benchmarks/benchmarks/analysis/rms.py
Expand Up @@ -96,13 +96,10 @@ def setup(self, n_atoms, step, weights):
self.u = MDAnalysis.Universe(PSF, DCD)
self.ag = self.u.atoms[:n_atoms]
self.RMSF_inst = rms.RMSF(atomgroup=self.ag,
start=None,
stop=None,
step=step,
weights=weights)

def time_RMSF(self, n_atoms, step, weights):
"""Benchmark RMSF.run() method, which parses
over the entire trajectory.
"""
self.RMSF_inst.run()
self.RMSF_inst.run(step=step)
11 changes: 10 additions & 1 deletion package/CHANGELOG
Expand Up @@ -78,7 +78,7 @@ Enhancements
analysis module (PR #1997, PR #2033)
* Added the analysis.data module for reference data used in analysis (PR #2033)
* Added analysis.dihedrals with Ramachandran class to analysis module (PR #1997)
* Added augment functionality to create relevant images of particles
* Added augment functionality to create relevant images of particles
* Most functions in `MDanalysis.lib.distances` previously only accepting
arrays of coordinates now also accept single coordinates as input (PR #2048,
Issues #1262 #1938)
Expand Down Expand Up @@ -109,6 +109,8 @@ Fixes
* Fixed SphericalLayer and SphericalZone selections with pbc=True. Previously
these shifted all atoms to inside the primary unit cell when calculating the
center of the reference group (Issue #1795)
* PCA analysis now uses start frame as reference frame rather than 0th frame
(PR #2055)

Changes
* TopologyAttrs are now statically typed (Issue #1876)
Expand All @@ -130,8 +132,15 @@ Changes
enforced to be of the form ``[lx, ly, lz, alpha, beta, gamma]`` as required
by the docs (Issue #2046, PR #2048)
* Added periodic keyword to select_atoms (#759)
* PCA.transform now requires that PCA.run() has already been called
otherwise a ValueError is raised (PR #2055)
* The quiet keyword has been removed and replaced with verbose (Issue #1975
PR #2055)

Deprecations
* start/stop/step are deprecated in the initialization of Analysis classes.
These parameters should instead be given to the run() method of the class.
(Issue #1463 #1979 PR #2055)
* Almost all "save()", "save_results()", "save_table()" methods in
analysis classes were deprecated and will be removed in 1.0 (Issue
#1972 and #1745)
Expand Down
88 changes: 54 additions & 34 deletions package/MDAnalysis/analysis/base.py
Expand Up @@ -32,11 +32,12 @@
from six.moves import range, zip
import inspect
import logging
import warnings

import numpy as np
from MDAnalysis import coordinates
from MDAnalysis.core.groups import AtomGroup
from MDAnalysis.lib.log import ProgressMeter, _set_verbose
from MDAnalysis.lib.log import ProgressMeter

logger = logging.getLogger(__name__)

Expand Down Expand Up @@ -84,30 +85,34 @@ def _conclude(self):
.. code-block:: python
na = NewAnalysis(u.select_atoms('name CA'), 35).run()
na = NewAnalysis(u.select_atoms('name CA'), 35).run(start=10, stop=20)
print(na.result)
"""

def __init__(self, trajectory, start=None,
stop=None, step=None, verbose=None, quiet=None):
def __init__(self, trajectory, verbose=False, **kwargs):
"""
Parameters
----------
trajectory : mda.Reader
A trajectory Reader
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
Turn on more logging and debugging, default ``False``
"""
self._verbose = _set_verbose(verbose, quiet, default=False)
self._quiet = not self._verbose
self._setup_frames(trajectory, start, stop, step)
self._trajectory = trajectory
self._verbose = verbose
# do deprecated kwargs
# remove in 1.0
deps = []
for arg in ['start', 'stop', 'step']:
if arg in kwargs and not kwargs[arg] is None:
deps.append(arg)
setattr(self, arg, kwargs[arg])
if deps:
warnings.warn('Setting the following kwargs should be '
'done in the run() method: {}'.format(
', '.join(deps)),
DeprecationWarning)

def _setup_frames(self, trajectory, start=None, stop=None, step=None):
"""
Expand All @@ -126,6 +131,11 @@ def _setup_frames(self, trajectory, start=None, stop=None, step=None):
number of frames to skip between each analysed frame
"""
self._trajectory = trajectory
# TODO: Remove once start/stop/step are deprecated from init
# See if these have been set as class attributes, and use that
start = getattr(self, 'start', start)
stop = getattr(self, 'stop', stop)
step = getattr(self, 'step', step)
start, stop, step = trajectory.check_slice_indices(start, stop, step)
self.start = start
self.stop = stop
Expand All @@ -135,21 +145,9 @@ def _setup_frames(self, trajectory, start=None, stop=None, step=None):
if interval == 0:
interval = 1

# ensure _verbose is set when __init__ wasn't called, this is to not
# break pre 0.16.0 API usage of AnalysisBase
if not hasattr(self, '_verbose'):
if hasattr(self, '_quiet'):
# Here, we are in the odd case where a children class defined
# self._quiet without going through AnalysisBase.__init__.
warnings.warn("The *_quiet* attribute of analyses is "
"deprecated (from 0.16)use *_verbose* instead.",
DeprecationWarning)
self._verbose = not self._quiet
else:
self._verbose = True
self._quiet = not self._verbose
verbose = getattr(self, '_verbose', False)
self._pm = ProgressMeter(self.n_frames if self.n_frames else 1,
interval=interval, verbose=self._verbose)
interval=interval, verbose=verbose)

def _single_frame(self):
"""Calculate data from a single frame of trajectory
Expand All @@ -169,8 +167,25 @@ def _conclude(self):
"""
pass

def run(self):
"""Perform the calculation"""
def run(self, start=None, stop=None, step=None, verbose=None):
"""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
"""
logger.info("Choosing frames to analyze")
# if verbose unchanged, use class default
verbose = getattr(self, '_verbose', False) if verbose is None else verbose

self._setup_frames(self._trajectory, start, stop, step)
logger.info("Starting preparation")
self._prepare()
for i, ts in enumerate(
Expand Down Expand Up @@ -241,10 +256,15 @@ def __init__(self, function, trajectory=None, *args, **kwargs):

self.function = function
self.args = args
base_kwargs, self.kwargs = _filter_baseanalysis_kwargs(self.function,
kwargs)

super(AnalysisFromFunction, self).__init__(trajectory, **base_kwargs)
# TODO: Remove in 1.0
my_kwargs = {}
for depped_arg in ['start', 'stop', 'step']:
if depped_arg in kwargs:
my_kwargs[depped_arg] = kwargs.pop(depped_arg)
self.kwargs = kwargs

super(AnalysisFromFunction, self).__init__(trajectory, **my_kwargs)

def _prepare(self):
self.results = []
Expand Down Expand Up @@ -272,7 +292,7 @@ def analysis_class(function):
>>> def RotationMatrix(mobile, ref):
>>> return mda.analysis.align.rotation_matrix(mobile, ref)[0]
>>> rot = RotationMatrix(u.trajectory, mobile, ref, step=2).run()
>>> rot = RotationMatrix(u.trajectory, mobile, ref).run(step=2)
>>> print(rot.results)
"""

Expand Down
8 changes: 0 additions & 8 deletions package/MDAnalysis/analysis/contacts.py
Expand Up @@ -388,14 +388,6 @@ def __init__(self, u, selection, refgroup, method="hard_cut", radius=4.5,
kwargs : dict, optional
dictionary of additional kwargs passed to `method`. Check
respective functions for reasonable values.
start : int, optional
First frame of trajectory to analyse, Default: None becomes 0.
stop : int, optional
Frame index to stop analysis. Default: None becomes
n_frames. Iteration stops *before* this frame number,
which means that the trajectory would be read until the end.
step : int, optional
Step between frames to analyse, Default: None becomes 1.
verbose : bool (optional)
Show detailed progress of the calculation if set to ``True``; the
default is ``False``.
Expand Down
4 changes: 2 additions & 2 deletions package/MDAnalysis/analysis/density.py
Expand Up @@ -552,7 +552,7 @@ def density_from_Universe(universe, delta=1.0, atomselection='name OH2',
start=None, stop=None, step=None,
metadata=None, padding=2.0, cutoff=0, soluteselection=None,
use_kdtree=True, update_selection=False,
verbose=None, interval=1, quiet=None,
verbose=False, interval=1, quiet=None,
parameters=None,
gridcenter=None, xdim=None, ydim=None, zdim=None):
"""Create a density grid from a :class:`MDAnalysis.Universe` object.
Expand Down Expand Up @@ -751,7 +751,7 @@ def current_coordinates():
h = grid.copy()

pm = ProgressMeter(u.trajectory.n_frames, interval=interval,
verbose=verbose, quiet=quiet,
verbose=verbose,
format="Histogramming %(n_atoms)6d atoms in frame "
"%(step)5d/%(numsteps)d [%(percentage)5.1f%%]\r")
start, stop, step = u.trajectory.check_slice_indices(start, stop, step)
Expand Down
44 changes: 23 additions & 21 deletions package/MDAnalysis/analysis/diffusionmap.py
Expand Up @@ -198,8 +198,7 @@ class DistanceMatrix(AnalysisBase):
"""
def __init__(self, u, select='all', metric=rmsd, cutoff=1E0-5,
weights=None, start=None, stop=None, step=None,
**kwargs):
weights=None, **kwargs):
"""
Parameters
----------
Expand All @@ -226,14 +225,6 @@ def __init__(self, u, select='all', metric=rmsd, cutoff=1E0-5,
Default: 1EO-5
weights : array, optional
Weights to be given to coordinates for metric calculation
start : int, optional
First frame of trajectory to analyse, Default: None becomes 0.
stop : int, optional
Frame index to stop analysis. Default: None becomes
n_frames. Iteration stops *before* this frame number,
which means that the trajectory would be read until the end.
step : int, optional
Step between frames to analyse, Default: None becomes 1.
verbose : bool (optional)
Show detailed progress of the calculation if set to ``True``; the
default is ``False``.
Expand All @@ -242,8 +233,7 @@ def __init__(self, u, select='all', metric=rmsd, cutoff=1E0-5,
traj = self._u.trajectory

# remember that this must be called before referencing self.n_frames
super(DistanceMatrix, self).__init__(self._u.trajectory,
start=start, stop=stop, step=step, **kwargs)
super(DistanceMatrix, self).__init__(self._u.trajectory, **kwargs)

self.atoms = self._u.select_atoms(select)
self._metric = metric
Expand Down Expand Up @@ -325,7 +315,7 @@ def __init__(self, u, epsilon=1, **kwargs):
**kwargs
Parameters to be passed for the initialization of a
:class:`DistanceMatrix`.
"""
"""
if isinstance(u, Universe):
self._dist_matrix = DistanceMatrix(u, **kwargs)
elif isinstance(u, DistanceMatrix):
Expand All @@ -334,20 +324,32 @@ def __init__(self, u, epsilon=1, **kwargs):
raise ValueError("U is not a Universe or DistanceMatrix and"
" so the DiffusionMap has no data to work with.")
self._epsilon = epsilon

def run(self, start=None, stop=None, step=None):
""" Create and decompose the diffusion matrix in preparation
for a diffusion map.
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
.. versionchanged:: 0.19.0
Added start/stop/step kwargs
"""
# run only if distance matrix not already calculated
if not self._dist_matrix._calculated:
self._dist_matrix.run(start=start, stop=stop, step=step)
# important for transform function and length of .run() method
self._n_frames = self._dist_matrix.n_frames
if self._n_frames > 5000:
warnings.warn("The distance matrix is very large, and can "
"be very slow to compute. Consider picking a larger "
"step size in distance matrix initialization.")


def run(self):
""" Create and decompose the diffusion matrix in preparation
for a diffusion map."""
# run only if distance matrix not already calculated
if not self._dist_matrix._calculated:
self._dist_matrix.run()
self._scaled_matrix = (self._dist_matrix.dist_matrix ** 2 /
self._epsilon)
# take negative exponent of scaled matrix to create Isotropic kernel
Expand Down

0 comments on commit 00cf3b7

Please sign in to comment.