forked from glotzerlab/freud
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtest_locality_Voronoi.py
417 lines (366 loc) · 13.6 KB
/
test_locality_Voronoi.py
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
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
import matplotlib
import matplotlib.pyplot as plt
import numpy as np
import numpy.testing as npt
import pytest
from util import sort_rounded_xyz_array
import freud
matplotlib.use("agg")
class TestVoronoi:
def test_random_2d(self):
# Test that voronoi tessellations of random systems have the same
# number of points and polytopes
L = 10 # Box length
N = 5000 # Number of particles
box, points = freud.data.make_random_system(L, N, is2D=True, seed=100)
vor = freud.locality.Voronoi()
vor.compute((box, points))
# Verify the polytopes and volumes
npt.assert_equal(len(vor.polytopes), len(points))
npt.assert_equal(len(vor.volumes), len(points))
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
# Ensure every point has neighbors
assert np.all(vor.nlist.neighbor_counts > 0)
# Ensure neighbor list is symmetric
ijs = set(zip(vor.nlist.query_point_indices, vor.nlist.point_indices))
jis = set(zip(vor.nlist.point_indices, vor.nlist.query_point_indices))
# Every (i, j) pair should have a corresponding (j, i) pair
assert all((j, i) in jis for (i, j) in ijs)
# The number of vertices in each polygon should be equal to
# the number of neighbors (only valid in 2D).
npt.assert_equal([len(p) for p in vor.polytopes], vor.nlist.neighbor_counts)
def test_random_3d(self):
# Test that voronoi tessellations of random systems have the same
# number of points and polytopes
L = 10 # Box length
N = 5000 # Number of particles
box, points = freud.data.make_random_system(L, N, is2D=False, seed=100)
vor = freud.locality.Voronoi()
vor.compute((box, points))
# Verify the polytopes and volumes
npt.assert_equal(len(vor.polytopes), len(points))
npt.assert_equal(len(vor.volumes), len(points))
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
# Ensure every point has neighbors
assert np.all(vor.nlist.neighbor_counts > 0)
# Ensure neighbor list is symmetric
ijs = set(zip(vor.nlist.query_point_indices, vor.nlist.point_indices))
jis = set(zip(vor.nlist.point_indices, vor.nlist.query_point_indices))
# Every (i, j) pair should have a corresponding (j, i) pair
assert all((j, i) in jis for (i, j) in ijs)
def test_voronoi_tess_2d(self):
# Test that the voronoi polytope works for a 2D system
L = 10 # Box length
box = freud.box.Box.square(L)
vor = freud.locality.Voronoi()
# Make a regular grid
points = np.array(
[
[0, 0, 0],
[0, 1, 0],
[0, 2, 0],
[1, 0, 0],
[1, 1, 0],
[1, 2, 0],
[2, 0, 0],
[2, 1, 0],
[2, 2, 0],
]
).astype(np.float64)
vor.compute((box, points))
center_polytope = sort_rounded_xyz_array(vor.polytopes[4])
expected_polytope = sort_rounded_xyz_array(
[[1.5, 1.5, 0], [0.5, 1.5, 0], [0.5, 0.5, 0], [1.5, 0.5, 0]]
)
npt.assert_almost_equal(center_polytope, expected_polytope)
# Verify the cell areas
npt.assert_almost_equal(vor.volumes[4], 1)
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
# Verify the neighbor list weights
npt.assert_almost_equal(
vor.nlist.weights[vor.nlist.query_point_indices == 4], 1
)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
# Double the points (still inside the box) and test again
points *= 2
vor.compute((box, points))
center_polytope = sort_rounded_xyz_array(vor.polytopes[4])
expected_polytope = sort_rounded_xyz_array(
[[3, 3, 0], [1, 3, 0], [1, 1, 0], [3, 1, 0]]
)
npt.assert_almost_equal(center_polytope, expected_polytope)
npt.assert_almost_equal(vor.volumes[4], 4)
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
npt.assert_almost_equal(
vor.nlist.weights[vor.nlist.query_point_indices == 4], 2
)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
def test_voronoi_tess_3d(self):
# Test that the voronoi polytope works for a 3D system
L = 10 # Box length
box = freud.box.Box.cube(L)
vor = freud.locality.Voronoi()
# Make a regular grid
points = np.array(
[
[0, 0, 0],
[0, 1, 0],
[0, 2, 0],
[1, 0, 0],
[1, 1, 0],
[1, 2, 0],
[2, 0, 0],
[2, 1, 0],
[2, 2, 0],
[0, 0, 1],
[0, 1, 1],
[0, 2, 1],
[1, 0, 1],
[1, 1, 1],
[1, 2, 1],
[2, 0, 1],
[2, 1, 1],
[2, 2, 1],
[0, 0, 2],
[0, 1, 2],
[0, 2, 2],
[1, 0, 2],
[1, 1, 2],
[1, 2, 2],
[2, 0, 2],
[2, 1, 2],
[2, 2, 2],
]
).astype(np.float64)
vor.compute((box, points))
center_polytope = sort_rounded_xyz_array(vor.polytopes[13])
expected_polytope = sort_rounded_xyz_array(
[
[1.5, 1.5, 1.5],
[1.5, 0.5, 1.5],
[1.5, 0.5, 0.5],
[1.5, 1.5, 0.5],
[0.5, 0.5, 0.5],
[0.5, 0.5, 1.5],
[0.5, 1.5, 0.5],
[0.5, 1.5, 1.5],
]
)
npt.assert_almost_equal(center_polytope, expected_polytope)
# Verify the cell volumes
npt.assert_almost_equal(vor.volumes[13], 1)
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
# Verify the neighbor list weights
npt.assert_almost_equal(
vor.nlist.weights[vor.nlist.query_point_indices == 13], 1
)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
# Double the points (still inside the box) and test again
points *= 2
vor.compute((box, points))
center_polytope = sort_rounded_xyz_array(vor.polytopes[13])
expected_polytope = sort_rounded_xyz_array(
[
[3, 3, 3],
[3, 1, 3],
[3, 1, 1],
[3, 3, 1],
[1, 1, 1],
[1, 1, 3],
[1, 3, 1],
[1, 3, 3],
]
)
npt.assert_almost_equal(center_polytope, expected_polytope)
npt.assert_almost_equal(vor.volumes[13], 8)
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
npt.assert_almost_equal(
vor.nlist.weights[vor.nlist.query_point_indices == 13], 4
)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
def test_voronoi_neighbors_wrapped(self):
# Test that voronoi neighbors in the first shell are correct for a
# wrapped 3D system, also tests multiple compute calls
n = 10
structure_neighbors = {
"sc": (freud.data.UnitCell.sc, 6),
"bcc": (freud.data.UnitCell.bcc, 14),
"fcc": (freud.data.UnitCell.fcc, 12),
}
vor = freud.locality.Voronoi()
for func, neighbors in structure_neighbors.values():
box, points = func().generate_system(n)
vor.compute((box, points))
nlist = vor.nlist
# Drop the tiny facets that come from numerical imprecision
nlist = nlist.filter(nlist.weights > 1e-5)
unique_indices, counts = np.unique(
nlist.query_point_indices, return_counts=True
)
# Every particle should have the specified number of neighbors
npt.assert_equal(counts, neighbors)
npt.assert_almost_equal(np.sum(vor.volumes), box.volume)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices]
- points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
def test_voronoi_weights_fcc(self):
# Test that voronoi neighbor weights are computed properly for 3D FCC
n = 3
box, points = freud.data.UnitCell.fcc().generate_system(n, scale=2)
vor = freud.locality.Voronoi()
vor.compute((box, points))
nlist = vor.nlist
# Drop the tiny facets that come from numerical imprecision
nlist = nlist.filter(nlist.weights > 1e-5)
# Every FCC particle should have 12 neighbors
npt.assert_equal(nlist.neighbor_counts, np.full(len(points), 12))
# Every facet area should be sqrt(2)/2
npt.assert_allclose(
nlist.weights, np.full(len(nlist.weights), 0.5 * np.sqrt(2)), atol=1e-5
)
# Every cell should have volume 2
vor.compute((box, points))
npt.assert_allclose(
vor.compute((box, points)).volumes,
np.full(len(vor.polytopes), 2.0),
atol=1e-5,
)
# Verify the neighbor distances
wrapped_distances = np.linalg.norm(
box.wrap(
points[vor.nlist.point_indices] - points[vor.nlist.query_point_indices]
),
axis=-1,
)
npt.assert_allclose(wrapped_distances, vor.nlist.distances)
def test_repr(self):
vor = freud.locality.Voronoi()
assert str(vor) == str(eval(repr(vor)))
def test_attributes(self):
# Test that the class attributes are protected
L = 10 # Box length
N = 40 # Number of particles
vor = freud.locality.Voronoi()
with pytest.raises(AttributeError):
vor.nlist
with pytest.raises(AttributeError):
vor.polytopes
with pytest.raises(AttributeError):
vor.volumes
box, points = freud.data.make_random_system(L, N, is2D=False)
vor.compute((box, points))
# Ensure attributes are accessible after calling compute
vor.nlist
vor.polytopes
vor.volumes
def test_repr_png(self):
L = 10 # Box length
box = freud.box.Box.square(L)
vor = freud.locality.Voronoi()
with pytest.raises(AttributeError):
vor.plot()
assert vor._repr_png_() is None
# Make a regular grid
points = np.array(
[
[0, 0, 0],
[0, 1, 0],
[0, 2, 0],
[1, 0, 0],
[1, 1, 0],
[1, 2, 0],
[2, 0, 0],
[2, 1, 0],
[2, 2, 0],
]
).astype(np.float32)
vor.compute((box, points))
vor._repr_png_()
L = 10 # Box length
box = freud.box.Box.cube(L)
vor = freud.locality.Voronoi()
# Make a regular grid
points = np.array(
[
[0, 0, 0],
[0, 1, 0],
[0, 2, 0],
[1, 0, 0],
[1, 1, 0],
[1, 2, 0],
[2, 0, 0],
[2, 1, 0],
[2, 2, 0],
[0, 0, 1],
[0, 1, 1],
[0, 2, 1],
[1, 0, 1],
[1, 1, 1],
[1, 2, 1],
[2, 0, 1],
[2, 1, 1],
[2, 2, 1],
[0, 0, 2],
[0, 1, 2],
[0, 2, 2],
[1, 0, 2],
[1, 1, 2],
[1, 2, 2],
[2, 0, 2],
[2, 1, 2],
[2, 2, 2],
]
).astype(np.float32)
vor.compute((box, points))
assert vor._repr_png_() is None
plt.close("all")