/usr/share/pyshared/ase/structure.py is in python-ase 3.6.0.2515-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
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 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 | """Atomic structure.
This mudule contains helper functions for setting up nanotubes and
graphene nanoribbons."""
import warnings
from math import sqrt
import numpy as np
from ase.atoms import Atoms, string2symbols
from ase.data import covalent_radii
from ase.utils import gcd
def nanotube(n, m, length=1, bond=1.42, symbol='C', verbose=False):
if n < m:
m, n = n, m
sign = -1
else:
sign = 1
nk = 6000
sq3 = sqrt(3.0)
a = sq3 * bond
l2 = n * n + m * m + n * m
l = sqrt(l2)
dt = a * l / np.pi
nd = gcd(n ,m)
if (n - m) % (3 * nd ) == 0:
ndr = 3 * nd
else:
ndr = nd
nr = (2 * m + n) / ndr
ns = -(2 * n + m) / ndr
nt2 = 3 * l2 / ndr / ndr
nt = np.floor(sqrt(nt2))
nn = 2 * l2 / ndr
ichk = 0
if nr == 0:
n60 = 1
else:
n60 = nr * 4
absn = abs(n60)
nnp = []
nnq = []
for i in range(-absn, absn + 1):
for j in range(-absn, absn + 1):
j2 = nr * j - ns * i
if j2 == 1:
j1 = m * i - n * j
if j1 > 0 and j1 < nn:
ichk += 1
nnp.append(i)
nnq.append(j)
if ichk == 0:
raise RuntimeError('not found p, q strange!!')
if ichk >= 2:
raise RuntimeError('more than 1 pair p, q strange!!')
nnnp = nnp[0]
nnnq = nnq[0]
if verbose:
print 'the symmetry vector is', nnnp, nnnq
lp = nnnp * nnnp + nnnq * nnnq + nnnp * nnnq
r = a * sqrt(lp)
c = a * l
t = sq3 * c / ndr
if 2 * nn > nk:
raise RuntimeError('parameter nk is too small!')
rs = c / (2.0 * np.pi)
if verbose:
print 'radius=', rs, t
q1 = np.arctan((sq3 * m) / (2 * n + m))
q2 = np.arctan((sq3 * nnnq) / (2 * nnnp + nnnq))
q3 = q1 - q2
q4 = 2.0 * np.pi / nn
q5 = bond * np.cos((np.pi / 6.0) - q1) / c * 2.0 * np.pi
h1 = abs(t) / abs(np.sin(q3))
h2 = bond * np.sin((np.pi / 6.0) - q1)
ii = 0
x, y, z = [], [], []
for i in range(nn):
x1, y1, z1 = 0, 0, 0
k = np.floor(i * abs(r) / h1)
x1 = rs * np.cos(i * q4)
y1 = rs * np.sin(i * q4)
z1 = (i * abs(r) - k * h1) * np.sin(q3)
kk2 = abs(np.floor((z1 + 0.0001) / t))
if z1 >= t - 0.0001:
z1 -= t * kk2
elif z1 < 0:
z1 += t * kk2
ii += 1
x.append(x1)
y.append(y1)
z.append(z1)
z3 = (i * abs(r) - k * h1) * np.sin(q3) - h2
ii += 1
if z3 >= 0 and z3 < t:
x2 = rs * np.cos(i * q4 + q5)
y2 = rs * np.sin(i * q4 + q5)
z2 = (i * abs(r) - k * h1) * np.sin(q3) - h2
x.append(x2)
y.append(y2)
z.append(z2)
else:
x2 = rs * np.cos(i * q4 + q5)
y2 = rs * np.sin(i * q4 + q5)
z2 = (i * abs(r) - (k + 1) * h1) * np.sin(q3) - h2
kk = abs(np.floor(z2 / t))
if z2 >= t - 0.0001:
z2 -= t * kk
elif z2 < 0:
z2 += t * kk
x.append(x2)
y.append(y2)
z.append(z2)
ntotal = 2 * nn
X = []
for i in range(ntotal):
X.append([x[i], y[i], sign * z[i]])
if length > 1:
xx = X[:]
for mnp in range(2, length + 1):
for i in range(len(xx)):
X.append(xx[i][:2] + [xx[i][2] + (mnp - 1) * t])
TransVec = t
NumAtom = ntotal * length
Diameter = rs * 2
ChiralAngle = np.arctan((sq3 * n) / (2 * m + n)) / (np.pi * 180)
cell = [Diameter * 2, Diameter * 2, length * t]
atoms = Atoms(symbol + str(NumAtom), positions=X, cell=cell,
pbc=[False, False, True])
atoms.center()
if verbose:
print 'translation vector =', TransVec
print 'diameter = ', Diameter
print 'chiral angle = ', ChiralAngle
return atoms
def graphene_nanoribbon(n, m, type='zigzag', saturated=False, C_H=1.09,
C_C=1.42, vacuum=2.5, magnetic=None, initial_mag=1.12,
sheet=False, main_element='C', saturate_element='H',
vacc=None):
"""Create a graphene nanoribbon.
Creates a graphene nanoribbon in the x-z plane, with the nanoribbon
running along the z axis.
Parameters:
n: The width of the nanoribbon
m: The length of the nanoribbon.
type ('zigzag'): The orientation of the ribbon. Must be either 'zigzag'
or 'armchair'.
saturated (Falsi): If true, hydrogen atoms are placed along the edge.
C_H: Carbon-hydrogen bond length. Default: 1.09 Angstrom
C_C: Carbon-carbon bond length. Default: 1.42 Angstrom.
vacuum: Amount of vacuum added to both sides. Default 2.5 Angstrom.
magnetic: Make the edges magnetic.
initial_mag: Magnitude of magnetic moment if magnetic=True.
sheet: If true, make an infinite sheet instead of a ribbon.
"""
#This function creates the coordinates for a graphene nanoribbon,
#n is width, m is length
if vacc is not None:
warnings.warn('Use vacuum=%f' % (0.5 * vacc))
vacuum = 0.5 * vacc
assert vacuum > 0
b = sqrt(3) * C_C / 4
arm_unit = Atoms(main_element+'4', pbc=(1,0,1),
cell = [4 * b, 2 * vacuum, 3 * C_C])
arm_unit.positions = [[0, 0, 0],
[b * 2, 0, C_C / 2.],
[b * 2, 0, 3 * C_C / 2.],
[0, 0, 2 * C_C]]
zz_unit = Atoms(main_element+'2', pbc=(1,0,1),
cell = [3 * C_C /2., 2 * vacuum, b * 4])
zz_unit.positions = [[0, 0, 0],
[C_C / 2., 0, b * 2]]
atoms = Atoms()
tol = 1e-4
if sheet:
vacuum2 = 0.0
else:
vacuum2 = vacuum
if type == 'zigzag':
edge_index0 = np.arange(m) * 2 + 1
edge_index1 = (n - 1) * m * 2 + np.arange(m) * 2
if magnetic:
mms = np.zeros(m * n * 2)
for i in edge_index0:
mms[i] = initial_mag
for i in edge_index1:
mms[i] = -initial_mag
for i in range(n):
layer = zz_unit.repeat((1, 1, m))
layer.positions[:, 0] -= 3 * C_C / 2 * i
if i % 2 == 1:
layer.positions[:, 2] += 2 * b
layer[-1].position[2] -= b * 4 * m
atoms += layer
if magnetic:
atoms.set_initial_magnetic_moments(mms)
if saturated:
H_atoms0 = Atoms(saturate_element + str(m))
H_atoms0.positions = atoms[edge_index0].positions
H_atoms0.positions[:, 0] += C_H
H_atoms1 = Atoms(saturate_element + str(m))
H_atoms1.positions = atoms[edge_index1].positions
H_atoms1.positions[:, 0] -= C_H
atoms += H_atoms0 + H_atoms1
atoms.cell = [n * 3 * C_C / 2 + 2 * vacuum2, 2 * vacuum, m * 4 * b]
elif type == 'armchair':
for i in range(n):
layer = arm_unit.repeat((1, 1, m))
layer.positions[:, 0] -= 4 * b * i
atoms += layer
atoms.cell = [b * 4 * n + 2 * vacuum2, 2 * vacuum, 3 * C_C * m]
atoms.center()
atoms.set_pbc([sheet, False, True])
return atoms
def molecule(name, data=None, **kwargs):
"""Create formula base on data. If data is None assume G2 set.
kwargs currently not used. """
if data is None:
from ase.data.g2 import data
if name not in data.keys():
raise NotImplementedError('%s not in data.' % (name))
args = data[name].copy()
# accept only the following Atoms constructor arguments
# XXX: should we accept all Atoms arguments?
for k in args.keys():
if k not in [
'symbols', 'positions', 'numbers',
'tags', 'masses',
'magmoms', 'charges',
'info',
]:
args.pop(k)
# kwargs overwrites data
args.update(kwargs)
return Atoms(**args)
def bulk(name, crystalstructure, a=None, c=None, covera=None,
orthorhombic=False, cubic=False):
"""Helper function for creating bulk systems.
name: str
Chemical symbol or symbols as in 'MgO' or 'NaCl'.
crystalstructure: str
Must be one of sc, fcc, bcc, hcp, diamond, zincblende or
rocksalt.
a: float
Lattice constant.
c: float
Lattice constant.
covera: float
c/a raitio used for hcp. Defaults to ideal ratio.
orthorhombic: bool
Construct orthorhombic unit cell instead of primitive cell
which is the default.
cubic: bool
Construct cubic unit cell.
"""
#warnings.warn('This function is deprecated. Use the ' +
# 'ase.lattice.bulk.bulk() function instead.')
if a is not None:
a = float(a)
if c is not None:
c = float(c)
if covera is not None and c is not None:
raise ValueError("Don't specify both c and c/a!")
if covera is None and c is None:
covera = sqrt(8.0 / 3.0)
if a is None:
a = estimate_lattice_constant(name, crystalstructure, covera)
if covera is None and c is not None:
covera = c / a
x = crystalstructure.lower()
if orthorhombic and x != 'sc':
return _orthorhombic_bulk(name, x, a, covera)
if cubic and x == 'bcc':
return _orthorhombic_bulk(name, x, a, covera)
if cubic and x != 'sc':
return _cubic_bulk(name, x, a)
if x == 'sc':
atoms = Atoms(name, cell=(a, a, a), pbc=True)
elif x == 'fcc':
b = a / 2
atoms = Atoms(name, cell=[(0, b, b), (b, 0, b), (b, b, 0)], pbc=True)
elif x == 'bcc':
b = a / 2
atoms = Atoms(name, cell=[(-b, b, b), (b, -b, b), (b, b, -b)],
pbc=True)
elif x == 'hcp':
atoms = Atoms(2 * name,
scaled_positions=[(0, 0, 0),
(1.0 / 3.0, 1.0 / 3.0, 0.5)],
cell=[(a, 0, 0),
(a / 2, a * sqrt(3) / 2, 0),
(0, 0, covera * a)],
pbc=True)
elif x == 'diamond':
atoms = bulk(2 * name, 'zincblende', a)
elif x == 'zincblende':
s1, s2 = string2symbols(name)
atoms = bulk(s1, 'fcc', a) + bulk(s2, 'fcc', a)
atoms.positions[1] += a / 4
elif x == 'rocksalt':
s1, s2 = string2symbols(name)
atoms = bulk(s1, 'fcc', a) + bulk(s2, 'fcc', a)
atoms.positions[1, 0] += a / 2
else:
raise ValueError('Unknown crystal structure: ' + crystalstructure)
return atoms
def estimate_lattice_constant(name, crystalstructure, covera):
atoms = bulk(name, crystalstructure, 1.0, covera)
v0 = atoms.get_volume()
v = 0.0
for Z in atoms.get_atomic_numbers():
r = covalent_radii[Z]
v += 4 * np.pi / 3 * r**3 * 1.5
return (v / v0)**(1.0 / 3)
def _orthorhombic_bulk(name, x, a, covera=None):
if x == 'fcc':
b = a / sqrt(2)
atoms = Atoms(2 * name, cell=(b, b, a), pbc=True,
scaled_positions=[(0, 0, 0), (0.5, 0.5, 0.5)])
elif x == 'bcc':
atoms = Atoms(2 * name, cell=(a, a, a), pbc=True,
scaled_positions=[(0, 0, 0), (0.5, 0.5, 0.5)])
elif x == 'hcp':
atoms = Atoms(4 * name,
cell=(a, a * sqrt(3), covera * a),
scaled_positions=[(0, 0, 0),
(0.5, 0.5, 0),
(0.5, 1.0 / 6.0, 0.5),
(0, 2.0 / 3.0, 0.5)],
pbc=True)
elif x == 'diamond':
atoms = _orthorhombic_bulk(2 * name, 'zincblende', a)
elif x == 'zincblende':
s1, s2 = string2symbols(name)
b = a / sqrt(2)
atoms = Atoms(2 * name, cell=(b, b, a), pbc=True,
scaled_positions=[(0, 0, 0), (0.5, 0, 0.25),
(0.5, 0.5, 0.5), (0, 0.5, 0.75)])
elif x == 'rocksalt':
s1, s2 = string2symbols(name)
b = a / sqrt(2)
atoms = Atoms(2 * name, cell=(b, b, a), pbc=True,
scaled_positions=[(0, 0, 0), (0.5, 0.5, 0),
(0.5, 0.5, 0.5), (0, 0, 0.5)])
else:
raise RuntimeError
return atoms
def _cubic_bulk(name, x, a):
if x == 'fcc':
atoms = Atoms(4 * name, cell=(a, a, a), pbc=True,
scaled_positions=[(0, 0, 0), (0, 0.5, 0.5),
(0.5, 0, 0.5), (0.5, 0.5, 0)])
elif x == 'diamond':
atoms = _cubic_bulk(2 * name, 'zincblende', a)
elif x == 'zincblende':
atoms = Atoms(4 * name, cell=(a, a, a), pbc=True,
scaled_positions=[(0, 0, 0), (0.25, 0.25, 0.25),
(0, 0.5, 0.5), (0.25, 0.75, 0.75),
(0.5, 0, 0.5), (0.75, 0.25, 0.75),
(0.5, 0.5, 0), (0.75, 0.75, 0.25)])
elif x == 'rocksalt':
atoms = Atoms(4 * name, cell=(a, a, a), pbc=True,
scaled_positions=[(0, 0, 0), (0.5, 0, 0),
(0, 0.5, 0.5), (0.5, 0.5, 0.5),
(0.5, 0, 0.5), (0, 0, 0.5),
(0.5, 0.5, 0), (0, 0.5, 0)])
else:
raise RuntimeError
return atoms
|