File: test_order.py

package info (click to toggle)
mdtraj 1.11.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 79,324 kB
  • sloc: python: 25,217; ansic: 6,266; cpp: 5,685; xml: 1,252; makefile: 192
file content (163 lines) | stat: -rwxr-xr-x 5,331 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
##############################################################################
# MDTraj: A Python Library for Loading, Saving, and Manipulating
#         Molecular Dynamics Trajectories.
# Copyright 2012-2017 Stanford University and the Authors
#
# Authors: Christoph Klein
# Contributors: Tim Moore
#
# MDTraj is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as
# published by the Free Software Foundation, either version 2.1
# of the License, or (at your option) any later version.
#
# This library is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public
# License along with MDTraj. If not, see <http://www.gnu.org/licenses/>.
##############################################################################

import numpy as np
import pytest

import mdtraj as md
from mdtraj.geometry import order
from mdtraj.testing import eq

"""The trajectories `2lines.pdb` and `3lines.pdb` contain several frames of,
respectively, 2 and 3 "residues" each consisting of two atoms in different
orientations relative to one another.

`2lines.pdb`
Frame 0: || in x         - S2 should be 1.0
Frame 1: || in y         - S2 should be 1.0
Frame 2: || in z         - S2 should be 1.0
Frame 3: |- in x/y       - S2 should be 0.25
Frame 4: |- in y/z       - S2 should be 0.25

`3lines.pdb`
Frame 0: ||| in x        - S2 should be 1.0
Frame 1: ||| in y        - S2 should be 1.0
Frame 2: ||| in z        - S2 should be 1.0
Frame 3: at right angles - S2 should be 0.0
"""


# TRAJ1 = md.load(get_fn('1line.pdb'))
# TRAJ2 = md.load(get_fn('2lines.pdb'))
# TRAJ2_LINE1 = TRAJ2.atom_slice((0, 1))
# TRAJ2_LINE2 = TRAJ2.atom_slice((2, 3))

# TRAJ3 = md.load(get_fn('3lines.pdb'))


@pytest.fixture()
def traj1(get_fn):
    return md.load(get_fn("1line.pdb"))


@pytest.fixture()
def traj2(get_fn):
    return md.load(get_fn("2lines.pdb"))


@pytest.fixture()
def traj3(get_fn):
    return md.load(get_fn("3lines.pdb"))


def test_nematic_order(traj2, traj3):
    assert eq(
        np.array([1.0, 1.0, 1.0, 0.25, 0.25]),
        order.compute_nematic_order(traj2),
    )
    assert eq(
        np.array([1.0, 1.0, 1.0, 0.0]),
        order.compute_nematic_order(traj3),
    )


def test_director(traj2):
    traj2_line1 = traj2.atom_slice((0, 1))
    traj2_line2 = traj2.atom_slice((2, 3))
    directors = order._compute_director(traj2_line1)
    assert eq(
        directors,
        np.array(
            [
                [1, 0, 0],  # Frame 0
                [0, 1, 0],  # Frame 1
                [0, 0, 1],  # Frame 2
                [1, 0, 0],  # Frame 3
                [1, 0, 0],  # Frame 4
            ],
        ),
    )
    directors = order._compute_director(traj2_line2)
    assert eq(
        directors,
        np.array(
            [
                [1, 0, 0],  # Frame 0
                [0, 1, 0],  # Frame 1
                [0, 0, 1],  # Frame 2
                [0, 1, 0],  # Frame 3
                [0, 0, 1],  # Frame 4
            ],
        ),
    )


def test_inertia(traj1, traj2, traj3):
    assert eq(
        order.compute_inertia_tensor(traj1),
        order._compute_inertia_tensor_slow(traj1),
    )
    assert eq(
        order.compute_inertia_tensor(traj2),
        order._compute_inertia_tensor_slow(traj2),
    )
    assert eq(
        order.compute_inertia_tensor(traj3),
        order._compute_inertia_tensor_slow(traj3),
    )


def test_nematic_order_args(traj2):
    with pytest.raises(ValueError):
        order.compute_nematic_order(traj2, indices="dog")
    with pytest.raises(ValueError):
        order.compute_nematic_order(traj2, indices="O")
    with pytest.raises(ValueError):
        order.compute_nematic_order(traj2, indices=1)

    # Input indices with wrong "shapes".
    with pytest.raises(ValueError):
        order.compute_nematic_order(traj2, indices=[[1, [2]], [2]])
    with pytest.raises(ValueError):
        order.compute_nematic_order(traj2, indices=[1, 2, 3])


def test_order_from_traj(get_fn):
    """Made a perfectly aligned monolayer, should have S2 = 1"""
    traj = md.load(get_fn("alkane-monolayer.pdb"))
    indices = [list(range(1900 + x, 1900 + x + 36)) for x in range(0, 64 * 36, 36)]
    s2 = md.compute_nematic_order(traj, indices=indices)
    np.testing.assert_allclose(1.0, s2)


def test_directors_angle_from_traj(get_fn):
    """All chains in example system aligned along z axis"""
    traj = md.load(get_fn("alkane-monolayer.pdb"))
    # only use the carbons for this calculation since they are excactly aligned
    # with the z-axis, and the hydrogens can throw things off just a bit
    indices = [list(range(1900 + x, 1900 + x + 30, 3)) for x in range(0, 64 * 36, 36)]
    directors = md.compute_directors(traj, indices=indices)
    for axis, target_angle in zip(np.eye(3), [90.0, 90.0, 0.0]):
        dotproducts = np.tensordot(directors, axis, axes=1)
        angles = np.degrees(np.arccos(dotproducts))
        angles = np.minimum(angles, 180 - angles)  # make vector point in positive z
        np.testing.assert_allclose(target_angle, angles, atol=1e-6)