root / ase / test / geometry.py @ 1
Historique | Voir | Annoter | Télécharger (4,18 ko)
1 | 1 | tkerber | "Test the ase.utils.geometry module"
|
---|---|---|---|
2 | 1 | tkerber | |
3 | 1 | tkerber | import numpy as np |
4 | 1 | tkerber | |
5 | 1 | tkerber | import ase |
6 | 1 | tkerber | from ase.lattice.spacegroup import crystal |
7 | 1 | tkerber | from ase.utils.geometry import get_layers, cut, stack |
8 | 1 | tkerber | |
9 | 1 | tkerber | np.set_printoptions(suppress=True)
|
10 | 1 | tkerber | |
11 | 1 | tkerber | al = crystal('Al', [(0,0,0)], spacegroup=225, cellpar=4.05) |
12 | 1 | tkerber | |
13 | 1 | tkerber | |
14 | 1 | tkerber | # Cut out slab of 5 Al(001) layers
|
15 | 1 | tkerber | al001 = cut(al, nlayers=5)
|
16 | 1 | tkerber | correct_pos = np.array([[ 0. , 0. , 0. ], |
17 | 1 | tkerber | [ 0. , 0.5, 0.2], |
18 | 1 | tkerber | [ 0.5, 0. , 0.2], |
19 | 1 | tkerber | [ 0.5, 0.5, 0. ], |
20 | 1 | tkerber | [ 0. , 0. , 0.4], |
21 | 1 | tkerber | [ 0. , 0.5, 0.6], |
22 | 1 | tkerber | [ 0.5, 0. , 0.6], |
23 | 1 | tkerber | [ 0.5, 0.5, 0.4], |
24 | 1 | tkerber | [ 0. , 0. , 0.8], |
25 | 1 | tkerber | [ 0.5, 0.5, 0.8]]) |
26 | 1 | tkerber | assert np.allclose(correct_pos, al001.get_scaled_positions())
|
27 | 1 | tkerber | |
28 | 1 | tkerber | # Check layers along 001
|
29 | 1 | tkerber | tags, levels = get_layers(al001, (0, 0, 1)) |
30 | 1 | tkerber | assert np.allclose(tags, [0, 1, 1, 0, 2, 3, 3, 2, 4, 4]) |
31 | 1 | tkerber | assert np.allclose(levels, [ 0., 2.025, 4.05, 6.075, 8.1]) |
32 | 1 | tkerber | |
33 | 1 | tkerber | # Check layers along 101
|
34 | 1 | tkerber | tags, levels = get_layers(al001, (1, 0, 1)) |
35 | 1 | tkerber | assert np.allclose(tags, [0, 1, 5, 3, 2, 4, 8, 7, 6, 9]) |
36 | 1 | tkerber | assert np.allclose(levels, [0.000, 0.752, 1.504, 1.880, 2.256, |
37 | 1 | tkerber | 2.632, 3.008, 3.384, 4.136, 4.888], atol=0.001) |
38 | 1 | tkerber | |
39 | 1 | tkerber | # Check layers along 111
|
40 | 1 | tkerber | tags, levels = get_layers(al001, (1, 1, 1)) |
41 | 1 | tkerber | assert np.allclose(tags, [0, 2, 2, 4, 1, 5, 5, 6, 3, 7]) |
42 | 1 | tkerber | assert np.allclose(levels, [0.000, 1.102, 1.929, 2.205, |
43 | 1 | tkerber | 2.756, 3.031, 3.858, 4.960], atol=0.001) |
44 | 1 | tkerber | |
45 | 1 | tkerber | |
46 | 1 | tkerber | # Cut out slab of three Al(111) layers
|
47 | 1 | tkerber | al111 = cut(al, (1,-1,0), (0,1,-1), nlayers=3) |
48 | 1 | tkerber | correct_pos = np.array([[ 0. , 0.5 , 0. ], |
49 | 1 | tkerber | [ 0.5 , 0.5 , 0. ], |
50 | 1 | tkerber | [ 2/3. , 5/6. , 1/3. ], |
51 | 1 | tkerber | [ 0.5 , 0. , 0. ], |
52 | 1 | tkerber | [ 0. , 0. , 0. ], |
53 | 1 | tkerber | [ 1/6. , 1/3. , 1/3. ], |
54 | 1 | tkerber | [ 2/3. , 1/3. , 1/3. ], |
55 | 1 | tkerber | [ 1/3. , 1/6. , 2/3. ], |
56 | 1 | tkerber | [ 5/6. , 1/6. , 2/3. ], |
57 | 1 | tkerber | [ 5/6. , 2/3. , 2/3. ], |
58 | 1 | tkerber | [ 1/6. , 5/6. , 1/3. ], |
59 | 1 | tkerber | [ 1/3. , 2/3. , 2/3. ]]) |
60 | 1 | tkerber | assert np.allclose(correct_pos, al111.get_scaled_positions())
|
61 | 1 | tkerber | |
62 | 1 | tkerber | |
63 | 1 | tkerber | # Cut out cell including all corner and edge atoms (non-periodic structure)
|
64 | 1 | tkerber | al = cut(al, extend=1.1)
|
65 | 1 | tkerber | correct_pos = np.array([[ 0. , 0. , 0. ], |
66 | 1 | tkerber | [ 0. , 0.5, 0.5], |
67 | 1 | tkerber | [ 0.5, 0. , 0.5], |
68 | 1 | tkerber | [ 0.5, 0.5, 0. ], |
69 | 1 | tkerber | [ 0. , 0. , 0. ], |
70 | 1 | tkerber | [ 0. , 0.5, 0.5], |
71 | 1 | tkerber | [ 0. , 0. , 0. ], |
72 | 1 | tkerber | [ 0.5, 0. , 0.5], |
73 | 1 | tkerber | [ 0. , 0. , 0. ], |
74 | 1 | tkerber | [ 0. , 0. , 0. ], |
75 | 1 | tkerber | [ 0.5, 0.5, 0. ], |
76 | 1 | tkerber | [ 0. , 0. , 0. ], |
77 | 1 | tkerber | [ 0. , 0. , 0. ], |
78 | 1 | tkerber | [ 0. , 0. , 0. ]]) |
79 | 1 | tkerber | assert np.allclose(correct_pos, al.get_scaled_positions())
|
80 | 1 | tkerber | |
81 | 1 | tkerber | |
82 | 1 | tkerber | # Create an Ag(111)/Si(111) interface
|
83 | 1 | tkerber | ag = crystal(['Ag'], basis=[(0,0,0)], spacegroup=225, cellpar=4.09) |
84 | 1 | tkerber | si = crystal(['Si'], basis=[(0,0,0)], spacegroup=227, cellpar=5.43) |
85 | 1 | tkerber | |
86 | 1 | tkerber | ag111 = cut(ag, a=(4, -4, 0), b=(4, 4, -8), nlayers=5) |
87 | 1 | tkerber | si111 = cut(si, a=(3, -3, 0), b=(3, 3, -6), nlayers=5) |
88 | 1 | tkerber | |
89 | 1 | tkerber | interface = stack(ag111, si111) |
90 | 1 | tkerber | assert len(interface) == 1000 |
91 | 1 | tkerber | assert np.allclose(interface.positions[::100], |
92 | 1 | tkerber | [[ 28.56875 , -2.040625 , -26.528125 ], |
93 | 1 | tkerber | [ 20.40916667, 6.12479167, -22.44395833], |
94 | 1 | tkerber | [ 24.49916667, 12.25541667, -20.39458333], |
95 | 1 | tkerber | [ 24.49041667, -6.11895833, -14.28145833], |
96 | 1 | tkerber | [ 10.20895833, 10.20895833, -12.23791667], |
97 | 1 | tkerber | [ 24.49625 , 2.049375 , -14.275625 ], |
98 | 1 | tkerber | [ 16.33375 , -4.0725 , 0.00875 ], |
99 | 1 | tkerber | [ 31.30027778, 12.25444444, -17.67472222], |
100 | 1 | tkerber | [ 25.85861111, 14.97527778, -14.95388889], |
101 | 1 | tkerber | [ 23.13777778, 4.09194444, -1.34972222]]) |