Skip to content

Commit

Permalink
documented selection of backend for distance calculations
Browse files Browse the repository at this point in the history
- closes issue #530.
- improved general docs for MDAnalysis.lib
  • Loading branch information
orbeckst committed Nov 12, 2015
1 parent d399b7b commit 30cbc19
Show file tree
Hide file tree
Showing 6 changed files with 76 additions and 22 deletions.
2 changes: 1 addition & 1 deletion package/MDAnalysis/lib/NeighborSearch.py
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,7 @@
#

"""
Neighbor Search wrapper for MDAnalysis --- :mod: `MDAnalysis.lib.NeighborSearch`
Neighbor Search wrapper for MDAnalysis --- :mod:`MDAnalysis.lib.NeighborSearch`
===============================================================================
This module contains classes that allow neighbor searches directly with
Expand Down
37 changes: 32 additions & 5 deletions package/MDAnalysis/lib/distances.py
Original file line number Diff line number Diff line change
Expand Up @@ -16,12 +16,14 @@
#

"""Fast distance array computation --- :mod:`MDAnalysis.lib.distances`
====================================================================
===================================================================
Fast C-routines to calculate distance arrays from coordinate arrays.
Fast C-routines to calculate distance arrays from coordinate
arrays. Many of the functions also exist in parallel versions, that
typically provide higher performance than the serial code.
Functions
---------
Selection of acceleration ("backend")
-------------------------------------
All functions take the optional keyword *backend*, which determines
the type of acceleration. Currently, the following choices are
Expand All @@ -36,6 +38,11 @@
with OpenMP
========== ======================== ======================================
.. versionadded:: 0.13.0
Functions
---------
.. autofunction:: distance_array(reference, configuration [, box [, result [, backend]]])
.. autofunction:: self_distance_array(reference [, box [,result [, backend]]])
.. autofunction:: calc_bonds(atom1, atom2 [, box, [, result [, backend]]])
Expand Down Expand Up @@ -219,7 +226,10 @@ def distance_array(reference, configuration, box=None, result=None, backend="ser
between reference coordinates i and configuration coordinates j
.. Note:: This method is slower than it could be because internally we need to
make copies of the ref and conf arrays.
make copies of the ref and conf arrays.
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
ref = reference.copy('C')
conf = configuration.copy('C')
Expand Down Expand Up @@ -300,6 +310,9 @@ def self_distance_array(reference, box=None, result=None, backend="serial"):
.. Note:: This method is slower than it could be because internally we need to
make copies of the coordinate arrays.
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
ref = reference.copy('C')

Expand Down Expand Up @@ -362,6 +375,9 @@ def transform_RtoS(inputcoords, box, backend="serial"):
:Returns:
*outcoords*
An n x 3 array of fractional coordiantes
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
coords = inputcoords.copy('C')
numcoords = coords.shape[0]
Expand Down Expand Up @@ -410,6 +426,9 @@ def transform_StoR(inputcoords, box, backend="serial"):
:Returns:
*outcoords*
An n x 3 array of fracional coordiantes
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
coords = inputcoords.copy('C')
numcoords = coords.shape[0]
Expand Down Expand Up @@ -473,6 +492,8 @@ def calc_bonds(coords1, coords2, box=None, result=None, backend="serial"):
numpy array with the length between each pair in coords1 and coords2
.. versionadded:: 0.8
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
atom1 = coords1.copy('C')
atom2 = coords2.copy('C')
Expand Down Expand Up @@ -558,6 +579,8 @@ def calc_angles(coords1, coords2, coords3, box=None, result=None, backend="seria
.. versionadded:: 0.8
.. versionchanged:: 0.9.0
Added optional box argument to account for periodic boundaries in calculation
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
atom1 = coords1.copy('C')
atom2 = coords2.copy('C')
Expand Down Expand Up @@ -656,6 +679,8 @@ def calc_dihedrals(coords1, coords2, coords3, coords4, box=None, result=None,
Added optional box argument to account for periodic boundaries in calculation
.. versionchanged:: 0.11.0
Renamed from calc_torsions to calc_dihedrals
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
atom1 = coords1.copy('C')
atom2 = coords2.copy('C')
Expand Down Expand Up @@ -727,6 +752,8 @@ def apply_PBC(incoords, box, backend="serial"):
as defined by box
.. versionadded:: 0.8
.. versionchanged:: 0.13.0
Added *backend* keyword.
"""
coords = incoords.copy('C')

Expand Down
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@
.. automodule:: MDAnalysis.lib._distances_openmp
:members:
16 changes: 16 additions & 0 deletions package/doc/sphinx/source/documentation_pages/lib/distances.rst
Original file line number Diff line number Diff line change
@@ -1,2 +1,18 @@
.. automodule:: MDAnalysis.lib.distances
:members:


Low-level modules for :mod:`MDAnalysis.lib.distances`
=====================================================

:mod:`MDAnalysis.lib._distances` contains low level access to the
*serial* MDAnalysis Cython functions in `distances`. These have
little to no error checking done on inputs so should be used with
caution. Similarly, :mod:`MDAnalysis.lib._distances_openmp` contains
the OpenMP-enable functions.

.. toctree::
:maxdepth: 1

_distances
_distances_openmp

This file was deleted.

39 changes: 25 additions & 14 deletions package/doc/sphinx/source/documentation_pages/lib_modules.rst
Original file line number Diff line number Diff line change
@@ -1,17 +1,27 @@
.. _lib:

************************
MDAnalysis library
************************
*******************************************
Library functions --- :mod:`MDAnalysis.lib`
*******************************************

:mod:`MDAnalysis.lib.distances` contains many high performance
maths functions.
.. module:: MDAnalysis.lib
:synopsis: ``lib`` collects independent code for re-use in MDAnalysis

:mod:`MDAnalysis.lib._distances` contains low level access to the
MDAnalysis' Cython functions in `distances`. These have little to no
error checking done on inputs so should be used with
caution. :mod:`MDAnalysis.lib.parallel` contains parallelized
implementations for some these algorithms.
:mod:`MDAnalysis.lib` contains code that is independent of the
specific MDAnalysis framework, such as fast calculators for distances
or simple logging support. Modules do not depend on other code inside
MDAnalysis except in :mod:`MDAnalysis.lib` itself (and possibly in
:mod:`MDAnalysis.exceptions`) and thus can be easily imported
elsewhere.

Overview
--------

:mod:`MDAnalysis.lib.distances` contains many high performance maths
functions. Most of them have the keyword *backend* that allows one to
either select serial (single threaded) code (``backend="serial``) or
to use parallelized versions (e.g. ``backend="OpenMP"`` for OpenMP
parallelism).

:mod:`MDAnalysis.lib.transformations` contains a multitude of
matrix operations for manipulating coordinate data.
Expand All @@ -22,19 +32,20 @@ superposition by minimizing the RMSD.
:mod:`MDAnalysis.lib.util` contains various file, string
and mathematical utility functions.

:mod:`MDAnalysis.lib.NeighborSearch` contains classes to do neighbor searches
with MDAnalysis objects.
:mod:`MDAnalysis.lib.NeighborSearch` contains classes to do neighbor
searches with MDAnalysis objects.


List of modules
---------------

.. toctree::
:maxdepth: 1

./lib/distances
./lib/_distances
./lib/NeighborSearch
./lib/log
./lib/mdamath
./lib/parallel
./lib/transformations
./lib/qcprot
./lib/util

0 comments on commit 30cbc19

Please sign in to comment.