File: test_neighbor.py

package info (click to toggle)
python-ase 3.21.1-2
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 13,936 kB
  • sloc: python: 122,428; xml: 946; makefile: 111; javascript: 47
file content (213 lines) | stat: -rw-r--r-- 7,283 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
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
import numpy.random as random
import numpy as np
import pytest
from ase import Atoms
from ase.neighborlist import (NeighborList, PrimitiveNeighborList,
                              NewPrimitiveNeighborList)
from ase.build import bulk


def count(nl, atoms):
    c = np.zeros(len(atoms), int)
    R = atoms.get_positions()
    cell = atoms.get_cell()
    d = 0.0
    for a in range(len(atoms)):
        i, offsets = nl.get_neighbors(a)
        for j in i:
            c[j] += 1
        c[a] += len(i)
        d += (((R[i] + np.dot(offsets, cell) - R[a])**2).sum(1)**0.5).sum()
    return d, c


# scipy sparse uses matrix subclass internally
@pytest.mark.filterwarnings('ignore:the matrix subclass')
@pytest.mark.slow
def test_supercell():
    atoms = Atoms(numbers=range(10),
                  cell=[(0.2, 1.2, 1.4),
                        (1.4, 0.1, 1.6),
                        (1.3, 2.0, -0.1)])
    atoms.set_scaled_positions(3 * random.random((10, 3)) - 1)

    for sorted in [False, True]:
        for p1 in range(2):
            for p2 in range(2):
                for p3 in range(2):
                    # print(p1, p2, p3)
                    atoms.set_pbc((p1, p2, p3))
                    nl = NeighborList(atoms.numbers * 0.2 + 0.5,
                                      skin=0.0, sorted=sorted)
                    nl.update(atoms)
                    d, c = count(nl, atoms)
                    atoms2 = atoms.repeat((p1 + 1, p2 + 1, p3 + 1))
                    nl2 = NeighborList(atoms2.numbers * 0.2 + 0.5,
                                       skin=0.0, sorted=sorted)
                    nl2.update(atoms2)
                    d2, c2 = count(nl2, atoms2)
                    c2.shape = (-1, 10)
                    dd = d * (p1 + 1) * (p2 + 1) * (p3 + 1) - d2
                    assert abs(dd) < 1e-10
                    assert not (c2 - c).any()


def test_H2():
    h2 = Atoms('H2', positions=[(0, 0, 0), (0, 0, 1)])
    nl = NeighborList([0.5, 0.5], skin=0.1, sorted=True, self_interaction=False)
    nl2 = NeighborList([0.5, 0.5], skin=0.1, sorted=True,
                       self_interaction=False,
                       primitive=NewPrimitiveNeighborList)
    assert nl2.update(h2)
    assert nl.update(h2)
    assert not nl.update(h2)
    assert (nl.get_neighbors(0)[0] == [1]).all()
    m = np.zeros((2, 2))
    m[0, 1] = 1
    assert np.array_equal(nl.get_connectivity_matrix(sparse=False), m)
    assert np.array_equal(nl.get_connectivity_matrix(sparse=True).todense(), m)
    assert np.array_equal(nl.get_connectivity_matrix().todense(),
                          nl2.get_connectivity_matrix().todense())

    h2[1].z += 0.09
    assert not nl.update(h2)
    assert (nl.get_neighbors(0)[0] == [1]).all()

    h2[1].z += 0.09
    assert nl.update(h2)
    assert (nl.get_neighbors(0)[0] == []).all()
    assert nl.nupdates == 2


def test_H2_shape_and_type():
    h2 = Atoms('H2', positions=[(0, 0, 0), (0, 0, 1)])
    nl = NeighborList([0.1, 0.1], skin=0.1, bothways=True,
                      self_interaction=False)
    assert nl.update(h2)
    assert nl.get_neighbors(0)[1].shape == (0, 3)
    assert nl.get_neighbors(0)[1].dtype == int


def test_fcc():
    x = bulk('X', 'fcc', a=2**0.5)

    nl = NeighborList([0.5], skin=0.01, bothways=True, self_interaction=False)
    nl.update(x)
    assert len(nl.get_neighbors(0)[0]) == 12

    nl = NeighborList([0.5] * 27, skin=0.01, bothways=True,
                      self_interaction=False)
    nl.update(x * (3, 3, 3))
    for a in range(27):
        assert len(nl.get_neighbors(a)[0]) == 12
    assert not np.any(nl.get_neighbors(13)[1])

    c = 0.0058
    for NeighborListClass in [PrimitiveNeighborList, NewPrimitiveNeighborList]:
        nl = NeighborListClass([c, c],
                               skin=0.0,
                               sorted=True,
                               self_interaction=False,
                               use_scaled_positions=True)
        nl.update([True, True, True],
                  np.eye(3) * 7.56,
                  np.array([[0, 0, 0],
                            [0, 0, 0.99875]]))
        n0, d0 = nl.get_neighbors(0)
        n1, d1 = nl.get_neighbors(1)
        # != is xor
        assert (np.all(n0 == [0]) and np.all(d0 == [0, 0, 1])) != \
            (np.all(n1 == [1]) and np.all(d1 == [0, 0, -1]))


def test_empty_neighbor_list():
    # Test empty neighbor list
    nl = PrimitiveNeighborList([])
    nl.update([True, True, True],
              np.eye(3) * 7.56,
              np.zeros((0, 3)))


def test_hexagonal_cell_and_large_cutoff():
    # Test hexagonal cell and large cutoff
    pbc_c = np.array([True, True, True])
    cutoff_a = np.array([8.0, 8.0])
    cell_cv = np.array([[0., 3.37316113, 3.37316113],
                        [3.37316113, 0., 3.37316113],
                        [3.37316113, 3.37316113, 0.]])
    spos_ac = np.array([[0., 0., 0.],
                        [0.25, 0.25, 0.25]])

    nl = PrimitiveNeighborList(cutoff_a, skin=0.0, sorted=True,
                               use_scaled_positions=True)
    nl2 = NewPrimitiveNeighborList(cutoff_a, skin=0.0, sorted=True,
                                   use_scaled_positions=True)
    nl.update(pbc_c, cell_cv, spos_ac)
    nl2.update(pbc_c, cell_cv, spos_ac)

    a0, offsets0 = nl.get_neighbors(0)
    b0 = np.zeros_like(a0)
    d0 = np.dot(spos_ac[a0] + offsets0 - spos_ac[0], cell_cv)
    a1, offsets1 = nl.get_neighbors(1)
    d1 = np.dot(spos_ac[a1] + offsets1 - spos_ac[1], cell_cv)
    b1 = np.ones_like(a1)

    a = np.concatenate([a0, a1])
    b = np.concatenate([b0, b1])
    d = np.concatenate([d0, d1])
    _a = np.concatenate([a, b])
    _b = np.concatenate([b, a])
    a = _a
    b = _b
    d = np.concatenate([d, -d])

    a0, offsets0 = nl2.get_neighbors(0)
    d0 = np.dot(spos_ac[a0] + offsets0 - spos_ac[0], cell_cv)
    b0 = np.zeros_like(a0)
    a1, offsets1 = nl2.get_neighbors(1)
    d1 = np.dot(spos_ac[a1] + offsets1 - spos_ac[1], cell_cv)
    b1 = np.ones_like(a1)

    a2 = np.concatenate([a0, a1])
    b2 = np.concatenate([b0, b1])
    d2 = np.concatenate([d0, d1])
    _a2 = np.concatenate([a2, b2])
    _b2 = np.concatenate([b2, a2])
    a2 = _a2
    b2 = _b2
    d2 = np.concatenate([d2, -d2])

    i = np.argsort(d[:, 0] + d[:, 1] * 1e2 + d[:, 2] * 1e4 + a * 1e6)
    i2 = np.argsort(d2[:, 0] + d2[:, 1] * 1e2 + d2[:, 2] * 1e4 + a2 * 1e6)

    assert np.all(a[i] == a2[i2])
    assert np.all(b[i] == b2[i2])
    assert np.allclose(d[i], d2[i2])


def test_small_cell_and_large_cutoff():
    # See: https://gitlab.com/ase/ase/-/issues/441
    cutoff = 50

    atoms = bulk('Cu', 'fcc', a=3.6)
    atoms *= (2, 2, 2)
    atoms.set_pbc(False)
    radii = cutoff * np.ones(len(atoms.get_atomic_numbers()))

    neighborhood_new = NeighborList(
        radii, skin=0.0, self_interaction=False, bothways=True,
        primitive=NewPrimitiveNeighborList
    )
    neighborhood_old = NeighborList(
        radii, skin=0.0, self_interaction=False, bothways=True,
        primitive=PrimitiveNeighborList
    )

    neighborhood_new.update(atoms)
    neighborhood_old.update(atoms)

    n0, d0 = neighborhood_new.get_neighbors(0)
    n1, d1 = neighborhood_old.get_neighbors(0)

    assert np.all(n0 == n1)
    assert np.all(d0 == d1)