1# -*- Mode: python; tab-width: 4; indent-tabs-mode:nil; coding:utf-8 -*-
2# vim: tabstop=4 expandtab shiftwidth=4 softtabstop=4
3#
4# MDAnalysis --- https://www.mdanalysis.org
5# Copyright (c) 2006-2017 The MDAnalysis Development Team and contributors
6# (see the file AUTHORS for the full list of names)
7#
8# Released under the GNU Public Licence, v2 or any higher version
9#
10# Please cite your use of MDAnalysis in published work:
11#
12# R. J. Gowers, M. Linke, J. Barnoud, T. J. E. Reddy, M. N. Melo, S. L. Seyler,
13# D. L. Dotson, J. Domanski, S. Buchoux, I. M. Kenney, and O. Beckstein.
14# MDAnalysis: A Python package for the rapid analysis of molecular dynamics
15# simulations. In S. Benthall and S. Rostrup editors, Proceedings of the 15th
16# Python in Science Conference, pages 102-109, Austin, TX, 2016. SciPy.
17#
18# N. Michaud-Agrawal, E. J. Denning, T. B. Woolf, and O. Beckstein.
19# MDAnalysis: A Toolkit for the Analysis of Molecular Dynamics Simulations.
20# J. Comput. Chem. 32 (2011), 2319--2327, doi:10.1002/jcc.21787
21#
22
23"""
24:mod:`MDAnalysis` --- analysis of molecular simulations in python
25=================================================================
26
27MDAnalysis (https://www.mdanalysis.org) is a python toolkit to analyze
28molecular dynamics trajectories generated by CHARMM, NAMD, Amber,
29Gromacs, or LAMMPS.
30
31It allows one to read molecular dynamics trajectories and access the
32atomic coordinates through numpy arrays. This provides a flexible and
33relatively fast framework for complex analysis tasks. In addition,
34CHARMM-style atom selection commands are implemented. Trajectories can
35also be manipulated (for instance, fit to a reference structure) and
36written out. Time-critical code is written in C for speed.
37
38Help is also available through the mailinglist at
39http://groups.google.com/group/mdnalysis-discussion
40
41Please report bugs and feature requests through the issue tracker at
42http://issues.mdanalysis.org
43
44Citation
45--------
46
47When using MDAnalysis in published work, please cite
48
49    N. Michaud-Agrawal, E. J. Denning, T. B. Woolf, and
50    O. Beckstein. MDAnalysis: A Toolkit for the Analysis of Molecular Dynamics
51    Simulations. J. Comput. Chem. 32 (2011), 2319--2327, doi:`10.1002/jcc.21787`_
52    https://www.mdanalysis.org
53
54For citations of included algorithms and sub-modules please see the references_.
55
56.. _`10.1002/jcc.21787`: http://dx.doi.org/10.1002/jcc.21787
57.. _references: https://docs.mdanalysis.org/documentation_pages/references.html
58
59
60Getting started
61---------------
62
63Import the package::
64
65  >>> import MDAnalysis
66
67(note that not everything in MDAnalysis is imported right away; for
68additional functionality you might have to import sub-modules
69separately, e.g. for RMS fitting ``import MDAnalysis.analysis.align``.)
70
71Build a "universe" from a topology (PSF, PDB) and a trajectory (DCD, XTC/TRR);
72here we are assuming that PSF, DCD, etc contain file names. If you don't have
73trajectories at hand you can play with the ones that come with MDAnalysis for
74testing (see below under `Examples`_)::
75
76  >>> u = MDAnalysis.Universe(PSF, DCD)
77
78Select the C-alpha atoms and store them as a group of atoms::
79
80  >>> ca = u.select_atoms('name CA')
81  >>> len(ca)
82  214
83
84Calculate the centre of mass of the CA and of all atoms::
85
86  >>> ca.center_of_mass()
87  array([ 0.06873595, -0.04605918, -0.24643682])
88  >>> u.atoms.center_of_mass()
89  array([-0.01094035,  0.05727601, -0.12885778])
90
91Calculate the CA end-to-end distance (in angstroem)::
92  >>> import numpy as np
93  >>> coord = ca.positions
94  >>> v = coord[-1] - coord[0]   # last Ca minus first one
95  >>> np.sqrt(np.dot(v, v,))
96  10.938133
97
98Define a function eedist():
99  >>> def eedist(atoms):
100  ...     coord = atoms.positions
101  ...     v = coord[-1] - coord[0]
102  ...     return sqrt(dot(v, v,))
103  ...
104  >>> eedist(ca)
105  10.938133
106
107and analyze all timesteps *ts* of the trajectory::
108  >>> for ts in u.trajectory:
109  ...      print eedist(ca)
110  10.9381
111  10.8459
112  10.4141
113   9.72062
114  ....
115
116See Also
117--------
118:class:`MDAnalysis.core.universe.Universe` for details
119
120
121Examples
122--------
123
124MDAnalysis comes with a number of real trajectories for testing. You
125can also use them to explore the functionality and ensure that
126everything is working properly::
127
128  from MDAnalysis import *
129  from MDAnalysis.tests.datafiles import PSF,DCD, PDB,XTC
130  u_dims_adk = Universe(PSF,DCD)
131  u_eq_adk = Universe(PDB, XTC)
132
133The PSF and DCD file are a closed-form-to-open-form transition of
134Adenylate Kinase (from [Beckstein2009]_) and the PDB+XTC file are ten
135frames from a Gromacs simulation of AdK solvated in TIP4P water with
136the OPLS/AA force field.
137
138.. [Beckstein2009] O. Beckstein, E.J. Denning, J.R. Perilla and T.B. Woolf,
139   Zipping and Unzipping of Adenylate Kinase: Atomistic Insights into the
140   Ensemble of Open <--> Closed Transitions. J Mol Biol 394 (2009), 160--176,
141   doi:10.1016/j.jmb.2009.09.009
142
143"""
144from __future__ import absolute_import
145
146__all__ = ['Universe', 'as_Universe', 'Writer', 'fetch_mmtf',
147           'AtomGroup', 'ResidueGroup', 'SegmentGroup']
148
149import logging
150import warnings
151
152logger = logging.getLogger("MDAnalysis.__init__")
153
154from .version import __version__
155try:
156    from .authors import __authors__
157except ImportError:
158    logger.info('Could not find authors.py, __authors__ will be empty.')
159    __authors__ = []
160
161# Registry of Readers, Parsers and Writers known to MDAnalysis
162# Metaclass magic fills these as classes are declared.
163_READERS = {}
164_SINGLEFRAME_WRITERS = {}
165_MULTIFRAME_WRITERS = {}
166_PARSERS = {}
167_SELECTION_WRITERS = {}
168# Registry of TopologyAttributes
169_TOPOLOGY_ATTRS = {}
170
171# Storing anchor universes for unpickling groups
172import weakref
173_ANCHOR_UNIVERSES = weakref.WeakValueDictionary()
174del weakref
175
176# custom exceptions and warnings
177from .exceptions import (
178    SelectionError, FinishTimeException, NoDataError, ApplicationError,
179    SelectionWarning, MissingDataWarning, ConversionWarning, FileFormatWarning,
180    StreamWarning
181)
182
183from .lib import log
184from .lib.log import start_logging, stop_logging
185
186logging.getLogger("MDAnalysis").addHandler(log.NullHandler())
187del logging
188
189# only MDAnalysis DeprecationWarnings are loud by default
190warnings.filterwarnings(action='once', category=DeprecationWarning,
191                        module='MDAnalysis')
192
193
194from . import units
195
196# Bring some often used objects into the current namespace
197from .core.universe import Universe, as_Universe, Merge
198from .core.groups import AtomGroup, ResidueGroup, SegmentGroup
199from .coordinates.core import writer as Writer
200
201# After Universe import
202from .coordinates.MMTF import fetch_mmtf
203
204from .migration.ten2eleven import ten2eleven
205
206from .due import due, Doi, BibTeX
207
208due.cite(BibTeX((
209            "@inproceedings{gowers2016, "
210            "title={MDAnalysis: A Python package for the rapid analysis "
211            "of molecular dynamics simulations}, "
212            "author={R. J. Gowers and M. Linke and "
213            "J. Barnoud and T. J. E. Reddy and M. N. Melo "
214            "and S. L. Seyler and D. L. Dotson and J. Domanski, and "
215            "S. Buchoux and I. M. Kenney and O. Beckstein},"
216            "journal={Proceedings of the 15th Python in Science Conference}, "
217            "pages={102-109}, "
218            "year={2016}, "
219            "editor={S. Benthall and S. Rostrup}, "
220            "note={Austin, TX, SciPy.} "
221            "}"
222            )),
223         description="Molecular simulation analysis library",
224         path="MDAnalysis", cite_module=True)
225due.cite(Doi("10.1002/jcc.21787"),
226         description="Molecular simulation analysis library",
227         path="MDAnalysis", cite_module=True)
228
229del Doi, BibTeX
230