root / ase / structure.py @ 16
Historique | Voir | Annoter | Télécharger (12,22 ko)
1 | 1 | tkerber | """Atomic structure.
|
---|---|---|---|
2 | 1 | tkerber |
|
3 | 1 | tkerber | This mudule contains helper functions for setting up nanotubes,
|
4 | 1 | tkerber | graphene nanoribbons and simple bulk crystals."""
|
5 | 1 | tkerber | |
6 | 1 | tkerber | from math import sqrt |
7 | 1 | tkerber | |
8 | 1 | tkerber | import numpy as np |
9 | 1 | tkerber | |
10 | 1 | tkerber | from ase.atoms import Atoms, string2symbols |
11 | 1 | tkerber | from ase.data import covalent_radii |
12 | 1 | tkerber | |
13 | 1 | tkerber | |
14 | 1 | tkerber | def gcd(m,n): |
15 | 1 | tkerber | while n:
|
16 | 1 | tkerber | m,n=n,m%n |
17 | 1 | tkerber | return m
|
18 | 1 | tkerber | |
19 | 1 | tkerber | def nanotube(n, m, length=1, bond=1.42, symbol='C', verbose=False): |
20 | 1 | tkerber | if n < m:
|
21 | 1 | tkerber | m, n = n, m |
22 | 1 | tkerber | nk = 6000
|
23 | 1 | tkerber | sq3 = sqrt(3.0)
|
24 | 1 | tkerber | a = sq3 * bond |
25 | 1 | tkerber | l2 = n * n + m * m + n * m |
26 | 1 | tkerber | l = sqrt(l2) |
27 | 1 | tkerber | dt = a * l / np.pi |
28 | 1 | tkerber | |
29 | 1 | tkerber | nd = gcd(n ,m) |
30 | 1 | tkerber | if (n - m) % (3 * nd ) == 0: |
31 | 1 | tkerber | ndr = 3 * nd
|
32 | 1 | tkerber | else:
|
33 | 1 | tkerber | ndr = nd |
34 | 1 | tkerber | |
35 | 1 | tkerber | nr = (2 * m + n) / ndr
|
36 | 1 | tkerber | ns = -(2 * n + m) / ndr
|
37 | 1 | tkerber | nt2 = 3 * l2 / ndr / ndr
|
38 | 1 | tkerber | nt = np.floor(sqrt(nt2)) |
39 | 1 | tkerber | nn = 2 * l2 / ndr
|
40 | 1 | tkerber | |
41 | 1 | tkerber | ichk = 0
|
42 | 1 | tkerber | if nr == 0: |
43 | 1 | tkerber | n60 = 1
|
44 | 1 | tkerber | else:
|
45 | 1 | tkerber | n60 = nr * 4
|
46 | 1 | tkerber | |
47 | 1 | tkerber | absn = abs(n60)
|
48 | 1 | tkerber | nnp = [] |
49 | 1 | tkerber | nnq = [] |
50 | 1 | tkerber | for i in range(-absn, absn + 1): |
51 | 1 | tkerber | for j in range(-absn, absn + 1): |
52 | 1 | tkerber | j2 = nr * j - ns * i |
53 | 1 | tkerber | if j2 == 1: |
54 | 1 | tkerber | j1 = m * i - n * j |
55 | 1 | tkerber | if j1 > 0 and j1 < nn: |
56 | 1 | tkerber | ichk += 1
|
57 | 1 | tkerber | nnp.append(i) |
58 | 1 | tkerber | nnq.append(j) |
59 | 1 | tkerber | |
60 | 1 | tkerber | if ichk == 0: |
61 | 1 | tkerber | raise RuntimeError('not found p, q strange!!') |
62 | 1 | tkerber | if ichk >= 2: |
63 | 1 | tkerber | raise RuntimeError('more than 1 pair p, q strange!!') |
64 | 1 | tkerber | |
65 | 1 | tkerber | nnnp = nnp[0]
|
66 | 1 | tkerber | nnnq = nnq[0]
|
67 | 1 | tkerber | |
68 | 1 | tkerber | if verbose:
|
69 | 1 | tkerber | print 'the symmetry vector is', nnnp, nnnq |
70 | 1 | tkerber | |
71 | 1 | tkerber | lp = nnnp * nnnp + nnnq * nnnq + nnnp * nnnq |
72 | 1 | tkerber | r = a * sqrt(lp) |
73 | 1 | tkerber | c = a * l |
74 | 1 | tkerber | t = sq3 * c / ndr |
75 | 1 | tkerber | |
76 | 1 | tkerber | if 2 * nn > nk: |
77 | 1 | tkerber | raise RuntimeError('parameter nk is too small!') |
78 | 1 | tkerber | |
79 | 1 | tkerber | rs = c / (2.0 * np.pi)
|
80 | 1 | tkerber | |
81 | 1 | tkerber | if verbose:
|
82 | 1 | tkerber | print 'radius=', rs, t |
83 | 1 | tkerber | |
84 | 1 | tkerber | q1 = np.arctan((sq3 * m) / (2 * n + m))
|
85 | 1 | tkerber | q2 = np.arctan((sq3 * nnnq) / (2 * nnnp + nnnq))
|
86 | 1 | tkerber | q3 = q1 - q2 |
87 | 1 | tkerber | |
88 | 1 | tkerber | q4 = 2.0 * np.pi / nn
|
89 | 1 | tkerber | q5 = bond * np.cos((np.pi / 6.0) - q1) / c * 2.0 * np.pi |
90 | 1 | tkerber | |
91 | 1 | tkerber | h1 = abs(t) / abs(np.sin(q3)) |
92 | 1 | tkerber | h2 = bond * np.sin((np.pi / 6.0) - q1)
|
93 | 1 | tkerber | |
94 | 1 | tkerber | ii = 0
|
95 | 1 | tkerber | x, y, z = [], [], [] |
96 | 1 | tkerber | for i in range(nn): |
97 | 1 | tkerber | x1, y1, z1 = 0, 0, 0 |
98 | 1 | tkerber | |
99 | 1 | tkerber | k = np.floor(i * abs(r) / h1)
|
100 | 1 | tkerber | x1 = rs * np.cos(i * q4) |
101 | 1 | tkerber | y1 = rs * np.sin(i * q4) |
102 | 1 | tkerber | z1 = (i * abs(r) - k * h1) * np.sin(q3)
|
103 | 1 | tkerber | kk2 = abs(np.floor((z1 + 0.0001) / t)) |
104 | 1 | tkerber | if z1 >= t - 0.0001: |
105 | 1 | tkerber | z1 -= t * kk2 |
106 | 1 | tkerber | elif z1 < 0: |
107 | 1 | tkerber | z1 += t * kk2 |
108 | 1 | tkerber | ii += 1
|
109 | 1 | tkerber | |
110 | 1 | tkerber | x.append(x1) |
111 | 1 | tkerber | y.append(y1) |
112 | 1 | tkerber | z.append(z1) |
113 | 1 | tkerber | z3 = (i * abs(r) - k * h1) * np.sin(q3) - h2
|
114 | 1 | tkerber | ii += 1
|
115 | 1 | tkerber | |
116 | 1 | tkerber | if z3 >= 0 and z3 < t: |
117 | 1 | tkerber | x2 = rs * np.cos(i * q4 + q5) |
118 | 1 | tkerber | y2 = rs * np.sin(i * q4 + q5) |
119 | 1 | tkerber | z2 = (i * abs(r) - k * h1) * np.sin(q3) - h2
|
120 | 1 | tkerber | x.append(x2) |
121 | 1 | tkerber | y.append(y2) |
122 | 1 | tkerber | z.append(z2) |
123 | 1 | tkerber | else:
|
124 | 1 | tkerber | x2 = rs * np.cos(i * q4 + q5) |
125 | 1 | tkerber | y2 = rs * np.sin(i * q4 + q5) |
126 | 1 | tkerber | z2 = (i * abs(r) - (k + 1) * h1) * np.sin(q3) - h2 |
127 | 1 | tkerber | kk = abs(np.floor(z2 / t))
|
128 | 1 | tkerber | if z2 >= t - 0.0001: |
129 | 1 | tkerber | z2 -= t * kk |
130 | 1 | tkerber | elif z2 < 0: |
131 | 1 | tkerber | z2 += t * kk |
132 | 1 | tkerber | x.append(x2) |
133 | 1 | tkerber | y.append(y2) |
134 | 1 | tkerber | z.append(z2) |
135 | 1 | tkerber | |
136 | 1 | tkerber | ntotal = 2 * nn
|
137 | 1 | tkerber | X = [] |
138 | 1 | tkerber | for i in range(ntotal): |
139 | 1 | tkerber | X.append([x[i], y[i], z[i]]) |
140 | 1 | tkerber | |
141 | 1 | tkerber | if length > 1: |
142 | 1 | tkerber | xx = X[:] |
143 | 1 | tkerber | for mnp in range(2, length + 1): |
144 | 1 | tkerber | for i in range(len(xx)): |
145 | 1 | tkerber | X.append(xx[i][:2] + [xx[i][2] + (mnp - 1) * t]) |
146 | 1 | tkerber | |
147 | 1 | tkerber | TransVec = t |
148 | 1 | tkerber | NumAtom = ntotal * length |
149 | 1 | tkerber | Diameter = rs * 2
|
150 | 1 | tkerber | ChiralAngle = np.arctan((sq3 * n) / (2 * m + n)) / (np.pi * 180) |
151 | 1 | tkerber | |
152 | 1 | tkerber | cell = [Diameter * 2, Diameter * 2, length * t] |
153 | 1 | tkerber | atoms = Atoms(symbol + str(NumAtom), positions=X, cell=cell,
|
154 | 1 | tkerber | pbc=[False, False, True]) |
155 | 1 | tkerber | atoms.center() |
156 | 1 | tkerber | if verbose:
|
157 | 1 | tkerber | print 'translation vector =', TransVec |
158 | 1 | tkerber | print 'diameter = ', Diameter |
159 | 1 | tkerber | print 'chiral angle = ', ChiralAngle |
160 | 1 | tkerber | return atoms
|
161 | 1 | tkerber | |
162 | 1 | tkerber | def graphene_nanoribbon(n, m, type='zigzag', saturated=False, C_H=1.09, |
163 | 1 | tkerber | C_C=1.42, vacc=5., magnetic=None, initial_mag=1.12, |
164 | 1 | tkerber | sheet=False, main_element='C', satur_element='H'): |
165 | 1 | tkerber | """Create a graphene nanoribbon.
|
166 | 1 | tkerber |
|
167 | 1 | tkerber | Creates a graphene nanoribbon in the x-z plane, with the nanoribbon
|
168 | 1 | tkerber | running along the z axis.
|
169 | 1 | tkerber |
|
170 | 1 | tkerber | Parameters:
|
171 | 1 | tkerber |
|
172 | 1 | tkerber | n: The width of the nanoribbon
|
173 | 1 | tkerber |
|
174 | 1 | tkerber | m: The length of the nanoribbon.
|
175 | 1 | tkerber |
|
176 | 1 | tkerber | type ('zigzag'): The orientation of the ribbon. Must be either 'zigzag'
|
177 | 1 | tkerber | or 'armchair'.
|
178 | 1 | tkerber |
|
179 | 1 | tkerber | saturated (Falsi): If true, hydrogen atoms are placed along the edge.
|
180 | 1 | tkerber |
|
181 | 1 | tkerber | C_H: Carbon-hydrogen bond length. Default: 1.09 Angstrom
|
182 | 1 | tkerber |
|
183 | 1 | tkerber | C_C: Carbon-carbon bond length. Default: 1.42 Angstrom.
|
184 | 1 | tkerber |
|
185 | 1 | tkerber | vacc: Amount of vacuum. Default 5 Angstrom.
|
186 | 1 | tkerber |
|
187 | 1 | tkerber | magnetic: Make the edges magnetic.
|
188 | 1 | tkerber |
|
189 | 1 | tkerber | initial_mag: Magnitude of magnetic moment if magnetic=True.
|
190 | 1 | tkerber |
|
191 | 1 | tkerber | sheet: If true, make an infinite sheet instead of a ribbon.
|
192 | 1 | tkerber | """
|
193 | 1 | tkerber | #This function creates the coordinates for a graphene nanoribbon,
|
194 | 1 | tkerber | #n is width, m is length
|
195 | 1 | tkerber | |
196 | 1 | tkerber | b = sqrt(3) * C_C / 4 |
197 | 1 | tkerber | arm_unit = Atoms(main_element+'4', pbc=(1,0,1), |
198 | 1 | tkerber | cell = [4 * b, vacc, 3 * C_C]) |
199 | 1 | tkerber | arm_unit.positions = [[0, 0, 0], |
200 | 1 | tkerber | [b * 2, 0, C_C / 2.], |
201 | 1 | tkerber | [b * 2, 0, 3 * C_C / 2.], |
202 | 1 | tkerber | [0, 0, 2 * C_C]] |
203 | 1 | tkerber | zz_unit = Atoms(main_element+'2', pbc=(1,0,1), |
204 | 1 | tkerber | cell = [3 * C_C /2., vacc, b * 4]) |
205 | 1 | tkerber | zz_unit.positions = [[0, 0, 0], |
206 | 1 | tkerber | [C_C / 2., 0, b * 2]] |
207 | 1 | tkerber | atoms = Atoms() |
208 | 1 | tkerber | tol = 1e-4
|
209 | 1 | tkerber | if sheet:
|
210 | 1 | tkerber | vacc2 = 0.0
|
211 | 1 | tkerber | else:
|
212 | 1 | tkerber | vacc2 = vacc |
213 | 1 | tkerber | if type == 'zigzag': |
214 | 1 | tkerber | edge_index0 = np.arange(m) * 2 + 1 |
215 | 1 | tkerber | edge_index1 = (n - 1) * m * 2 + np.arange(m) * 2 |
216 | 1 | tkerber | if magnetic:
|
217 | 1 | tkerber | mms = np.zeros(m * n * 2)
|
218 | 1 | tkerber | for i in edge_index0: |
219 | 1 | tkerber | mms[i] = initial_mag |
220 | 1 | tkerber | for i in edge_index1: |
221 | 1 | tkerber | mms[i] = -initial_mag |
222 | 1 | tkerber | |
223 | 1 | tkerber | for i in range(n): |
224 | 1 | tkerber | layer = zz_unit.repeat((1, 1, m)) |
225 | 1 | tkerber | layer.positions[:, 0] -= 3 * C_C / 2 * i |
226 | 1 | tkerber | if i % 2 == 1: |
227 | 1 | tkerber | layer.positions[:, 2] += 2 * b |
228 | 1 | tkerber | layer[-1].position[2] -= b * 4 * m |
229 | 1 | tkerber | atoms += layer |
230 | 1 | tkerber | if magnetic:
|
231 | 1 | tkerber | atoms.set_initial_magnetic_moments(mms) |
232 | 1 | tkerber | if saturated:
|
233 | 1 | tkerber | H_atoms0 = Atoms(satur_element + str(m))
|
234 | 1 | tkerber | H_atoms0.positions = atoms[edge_index0].positions |
235 | 1 | tkerber | H_atoms0.positions[:, 0] += C_H
|
236 | 1 | tkerber | H_atoms1 = Atoms(satur_element + str(m))
|
237 | 1 | tkerber | H_atoms1.positions = atoms[edge_index1].positions |
238 | 1 | tkerber | H_atoms1.positions[:, 0] -= C_H
|
239 | 1 | tkerber | atoms += H_atoms0 + H_atoms1 |
240 | 1 | tkerber | atoms.cell = [n * 3 * C_C / 2 + vacc2, vacc, m * 4 * b] |
241 | 1 | tkerber | |
242 | 1 | tkerber | elif type == 'armchair': |
243 | 1 | tkerber | for i in range(n): |
244 | 1 | tkerber | layer = arm_unit.repeat((1, 1, m)) |
245 | 1 | tkerber | layer.positions[:, 0] -= 4 * b * i |
246 | 1 | tkerber | atoms += layer |
247 | 1 | tkerber | atoms.cell = [b * 4 * n + vacc2, vacc, 3 * C_C * m] |
248 | 1 | tkerber | atoms.center() |
249 | 1 | tkerber | atoms.set_pbc([sheet, False, True]) |
250 | 1 | tkerber | return atoms
|
251 | 1 | tkerber | |
252 | 1 | tkerber | |
253 | 1 | tkerber | def bulk(name, crystalstructure, a=None, c=None, covera=None, |
254 | 1 | tkerber | orthorhombic=False, cubic=False): |
255 | 1 | tkerber | """Helper function for creating bulk systems.
|
256 | 1 | tkerber |
|
257 | 1 | tkerber | name: str
|
258 | 1 | tkerber | Chemical symbol or symbols as in 'MgO' or 'NaCl'.
|
259 | 1 | tkerber | crystalstructure: str
|
260 | 1 | tkerber | Must be one of sc, fcc, bcc, hcp, diamond, zincblende or
|
261 | 1 | tkerber | rocksalt.
|
262 | 1 | tkerber | a: float
|
263 | 1 | tkerber | Lattice constant.
|
264 | 1 | tkerber | c: float
|
265 | 1 | tkerber | Lattice constant.
|
266 | 1 | tkerber | covera: float
|
267 | 1 | tkerber | c/a raitio used for hcp. Defaults to ideal ratio.
|
268 | 1 | tkerber | orthorhombic: bool
|
269 | 1 | tkerber | Construct orthorhombic unit cell instead of primitive cell
|
270 | 1 | tkerber | which is the default.
|
271 | 1 | tkerber | cubic: bool
|
272 | 1 | tkerber | Construct cubic unit cell.
|
273 | 1 | tkerber | """
|
274 | 1 | tkerber | |
275 | 1 | tkerber | if covera is not None and c is not None: |
276 | 1 | tkerber | raise ValueError("Don't specify both c and c/a!") |
277 | 1 | tkerber | |
278 | 1 | tkerber | if covera is None and c is None: |
279 | 1 | tkerber | covera = sqrt(8.0 / 3.0) |
280 | 1 | tkerber | |
281 | 1 | tkerber | if a is None: |
282 | 1 | tkerber | a = estimate_lattice_constant(name, crystalstructure, covera) |
283 | 1 | tkerber | |
284 | 1 | tkerber | if covera is None and c is not None: |
285 | 1 | tkerber | covera = c / a |
286 | 1 | tkerber | |
287 | 1 | tkerber | x = crystalstructure.lower() |
288 | 1 | tkerber | |
289 | 1 | tkerber | if orthorhombic and x != 'sc': |
290 | 1 | tkerber | return _orthorhombic_bulk(name, x, a, covera)
|
291 | 1 | tkerber | |
292 | 1 | tkerber | if cubic and x == 'bcc': |
293 | 1 | tkerber | return _orthorhombic_bulk(name, x, a, covera)
|
294 | 1 | tkerber | |
295 | 1 | tkerber | if cubic and x != 'sc': |
296 | 1 | tkerber | return _cubic_bulk(name, x, a)
|
297 | 1 | tkerber | |
298 | 1 | tkerber | if x == 'sc': |
299 | 1 | tkerber | atoms = Atoms(name, cell=(a, a, a), pbc=True)
|
300 | 1 | tkerber | elif x == 'fcc': |
301 | 1 | tkerber | b = a / 2
|
302 | 1 | tkerber | atoms = Atoms(name, cell=[(0, b, b), (b, 0, b), (b, b, 0)], pbc=True) |
303 | 1 | tkerber | elif x == 'bcc': |
304 | 1 | tkerber | b = a / 2
|
305 | 1 | tkerber | atoms = Atoms(name, cell=[(-b, b, b), (b, -b, b), (b, b, -b)], |
306 | 1 | tkerber | pbc=True)
|
307 | 1 | tkerber | elif x == 'hcp': |
308 | 1 | tkerber | atoms = Atoms(2 * name,
|
309 | 1 | tkerber | scaled_positions=[(0, 0, 0), |
310 | 1 | tkerber | (1.0 / 3.0, 1.0 / 3.0, 0.5)], |
311 | 1 | tkerber | cell=[(a, 0, 0), |
312 | 1 | tkerber | (a / 2, a * sqrt(3) / 2, 0), |
313 | 1 | tkerber | (0, 0, covera * a)], |
314 | 1 | tkerber | pbc=True)
|
315 | 1 | tkerber | elif x == 'diamond': |
316 | 1 | tkerber | atoms = bulk(2 * name, 'zincblende', a) |
317 | 1 | tkerber | elif x == 'zincblende': |
318 | 1 | tkerber | s1, s2 = string2symbols(name) |
319 | 1 | tkerber | atoms = bulk(s1, 'fcc', a) + bulk(s2, 'fcc', a) |
320 | 1 | tkerber | atoms.positions[1] += a / 4 |
321 | 1 | tkerber | elif x == 'rocksalt': |
322 | 1 | tkerber | s1, s2 = string2symbols(name) |
323 | 1 | tkerber | atoms = bulk(s1, 'fcc', a) + bulk(s2, 'fcc', a) |
324 | 1 | tkerber | atoms.positions[1, 0] += a / 2 |
325 | 1 | tkerber | else:
|
326 | 1 | tkerber | raise ValueError('Unknown crystal structure: ' + crystalstructure) |
327 | 1 | tkerber | |
328 | 1 | tkerber | return atoms
|
329 | 1 | tkerber | |
330 | 1 | tkerber | def estimate_lattice_constant(name, crystalstructure, covera): |
331 | 1 | tkerber | atoms = bulk(name, crystalstructure, 1.0, covera)
|
332 | 1 | tkerber | v0 = atoms.get_volume() |
333 | 1 | tkerber | v = 0.0
|
334 | 1 | tkerber | for Z in atoms.get_atomic_numbers(): |
335 | 1 | tkerber | r = covalent_radii[Z] |
336 | 1 | tkerber | v += 4 * np.pi / 3 * r**3 * 1.5 |
337 | 1 | tkerber | return (v / v0)**(1.0 / 3) |
338 | 1 | tkerber | |
339 | 1 | tkerber | def _orthorhombic_bulk(name, x, a, covera=None): |
340 | 1 | tkerber | if x == 'fcc': |
341 | 1 | tkerber | b = a / sqrt(2)
|
342 | 1 | tkerber | atoms = Atoms(2 * name, cell=(b, b, a), pbc=True, |
343 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0.5, 0.5, 0.5)]) |
344 | 1 | tkerber | elif x == 'bcc': |
345 | 1 | tkerber | atoms = Atoms(2 * name, cell=(a, a, a), pbc=True, |
346 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0.5, 0.5, 0.5)]) |
347 | 1 | tkerber | elif x == 'hcp': |
348 | 1 | tkerber | atoms = Atoms(4 * name,
|
349 | 1 | tkerber | cell=(a, a * sqrt(3), covera * a),
|
350 | 1 | tkerber | scaled_positions=[(0, 0, 0), |
351 | 1 | tkerber | (0.5, 0.5, 0), |
352 | 1 | tkerber | (0.5, 1.0 / 6.0, 0.5), |
353 | 1 | tkerber | (0, 2.0 / 3.0, 0.5)], |
354 | 1 | tkerber | pbc=True)
|
355 | 1 | tkerber | elif x == 'diamond': |
356 | 1 | tkerber | atoms = _orthorhombic_bulk(2 * name, 'zincblende', a) |
357 | 1 | tkerber | elif x == 'zincblende': |
358 | 1 | tkerber | s1, s2 = string2symbols(name) |
359 | 1 | tkerber | b = a / sqrt(2)
|
360 | 1 | tkerber | atoms = Atoms(2 * name, cell=(b, b, a), pbc=True, |
361 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0.5, 0, 0.25), |
362 | 1 | tkerber | (0.5, 0.5, 0.5), (0, 0.5, 0.75)]) |
363 | 1 | tkerber | elif x == 'rocksalt': |
364 | 1 | tkerber | s1, s2 = string2symbols(name) |
365 | 1 | tkerber | b = a / sqrt(2)
|
366 | 1 | tkerber | atoms = Atoms(2 * name, cell=(b, b, a), pbc=True, |
367 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0.5, 0.5, 0), |
368 | 1 | tkerber | (0.5, 0.5, 0.5), (0, 0, 0.5)]) |
369 | 1 | tkerber | else:
|
370 | 1 | tkerber | raise RuntimeError |
371 | 1 | tkerber | |
372 | 1 | tkerber | return atoms
|
373 | 1 | tkerber | |
374 | 1 | tkerber | def _cubic_bulk(name, x, a): |
375 | 1 | tkerber | if x == 'fcc': |
376 | 1 | tkerber | atoms = Atoms(4 * name, cell=(a, a, a), pbc=True, |
377 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0, 0.5, 0.5), |
378 | 1 | tkerber | (0.5, 0, 0.5), (0.5, 0.5, 0)]) |
379 | 1 | tkerber | elif x == 'diamond': |
380 | 1 | tkerber | atoms = _cubic_bulk(2 * name, 'zincblende', a) |
381 | 1 | tkerber | elif x == 'zincblende': |
382 | 1 | tkerber | atoms = Atoms(4 * name, cell=(a, a, a), pbc=True, |
383 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0.25, 0.25, 0.25), |
384 | 1 | tkerber | (0, 0.5, 0.5), (0.25, 0.75, 0.75), |
385 | 1 | tkerber | (0.5, 0, 0.5), (0.75, 0.25, 0.75), |
386 | 1 | tkerber | (0.5, 0.5, 0), (0.75, 0.75, 0.25)]) |
387 | 1 | tkerber | elif x == 'rocksalt': |
388 | 1 | tkerber | atoms = Atoms(4 * name, cell=(a, a, a), pbc=True, |
389 | 1 | tkerber | scaled_positions=[(0, 0, 0), (0.5, 0, 0), |
390 | 1 | tkerber | (0, 0.5, 0.5), (0.5, 0.5, 0.5), |
391 | 1 | tkerber | (0.5, 0, 0.5), (0, 0, 0.5), |
392 | 1 | tkerber | (0.5, 0.5, 0), (0, 0.5, 0)]) |
393 | 1 | tkerber | else:
|
394 | 1 | tkerber | raise RuntimeError |
395 | 1 | tkerber | |
396 | 1 | tkerber | return atoms |