Statistiques
| Révision :

root / ase / utils / geometry.py @ 18

Historique | Voir | Annoter | Télécharger (12,3 ko)

1 1 tkerber
# Copyright (C) 2010, Jesper Friis
2 1 tkerber
# (see accompanying license files for details).
3 1 tkerber
4 1 tkerber
"""Utility tools for convenient creation of slabs and interfaces of
5 1 tkerber
different orientations."""
6 1 tkerber
7 1 tkerber
import numpy as np
8 1 tkerber
9 1 tkerber
10 1 tkerber
11 1 tkerber
def gcd(seq):
12 1 tkerber
    """Returns greatest common divisor of integers in *seq*."""
13 1 tkerber
    def _gcd(m, n):
14 1 tkerber
        while n:
15 1 tkerber
            m, n = n, m%n
16 1 tkerber
        return m
17 1 tkerber
    return reduce(_gcd, seq)
18 1 tkerber
19 1 tkerber
20 1 tkerber
21 1 tkerber
22 1 tkerber
def get_layers(atoms, miller, tolerance=0.001):
23 1 tkerber
    """Returns two arrays describing which layer each atom belongs
24 1 tkerber
    to and the distance between the layers and origo.
25 1 tkerber

26 1 tkerber
    Parameters:
27 1 tkerber

28 1 tkerber
    miller: 3 integers
29 1 tkerber
        The Miller indices of the planes. Actually, any direction
30 1 tkerber
        in reciprocal space works, so if a and b are two float
31 1 tkerber
        vectors spanning an atomic plane, you can get all layers
32 1 tkerber
        parallel to this with miller=np.cross(a,b).
33 1 tkerber
    tolerance: float
34 1 tkerber
        The maximum distance in Angstrom along the plane normal for
35 1 tkerber
        counting two atoms as belonging to the same plane.
36 1 tkerber

37 1 tkerber
    Returns:
38 1 tkerber

39 1 tkerber
    tags: array of integres
40 1 tkerber
        Array of layer indices for each atom.
41 1 tkerber
    levels: array of floats
42 1 tkerber
        Array of distances in Angstrom from each layer to origo.
43 1 tkerber

44 1 tkerber
    Example:
45 1 tkerber

46 1 tkerber
    >>> import numpy as np
47 1 tkerber
    >>> from ase.lattice.spacegroup import crystal
48 1 tkerber
    >>> atoms = crystal('Al', [(0,0,0)], spacegroup=225, cellpar=4.05)
49 1 tkerber
    >>> np.round(atoms.positions, decimals=5)
50 1 tkerber
    array([[ 0.   ,  0.   ,  0.   ],
51 1 tkerber
           [ 0.   ,  2.025,  2.025],
52 1 tkerber
           [ 2.025,  0.   ,  2.025],
53 1 tkerber
           [ 2.025,  2.025,  0.   ]])
54 1 tkerber
    >>> get_layers(atoms, (0,0,1))
55 1 tkerber
    (array([0, 1, 1, 0]), array([ 0.   ,  2.025]))
56 1 tkerber
    """
57 1 tkerber
    miller = np.asarray(miller)
58 1 tkerber
59 1 tkerber
    metric = np.dot(atoms.cell, atoms.cell.T)
60 1 tkerber
    c = np.linalg.solve(metric.T, miller.T).T
61 1 tkerber
    miller_norm = np.sqrt(np.dot(c, miller))
62 1 tkerber
    d = np.dot(atoms.get_scaled_positions(), miller)/miller_norm
63 1 tkerber
64 1 tkerber
    keys = np.argsort(d)
65 1 tkerber
    ikeys = np.argsort(keys)
66 1 tkerber
    mask = np.concatenate(([True], np.diff(d[keys]) > tolerance))
67 1 tkerber
    tags = np.cumsum(mask)[ikeys]
68 1 tkerber
    if tags.min() == 1:
69 1 tkerber
        tags -= 1
70 1 tkerber
71 1 tkerber
    levels = d[keys][mask]
72 1 tkerber
    return tags, levels
73 1 tkerber
74 1 tkerber
75 1 tkerber
76 1 tkerber
77 1 tkerber
def cut(atoms, a=(1,0,0), b=(0,1,0), c=(0,0,1), origo=(0,0,0),
78 1 tkerber
        nlayers=None, extend=1.0, tolerance=0.001):
79 1 tkerber
    """Cuts out a cell defined by *a*, *b*, *c* and *origo* from a
80 1 tkerber
    sufficiently repeated copy of *atoms*.
81 1 tkerber

82 1 tkerber
    Typically, this function is used to create slabs of different
83 1 tkerber
    sizes and orientations. The vectors *a*, *b* and *c* are in scaled
84 1 tkerber
    coordinates and defines the returned cell and should normally be
85 1 tkerber
    integer-valued in order to end up with a periodic
86 1 tkerber
    structure. However, for systems with sub-translations, like fcc,
87 1 tkerber
    integer multiples of 1/2 or 1/3 might also make sence for some
88 1 tkerber
    directions (and will be treated correctly).
89 1 tkerber

90 1 tkerber
    Parameters:
91 1 tkerber

92 1 tkerber
    atoms: Atoms instance
93 1 tkerber
        This should correspond to a repeatable unit cell.
94 1 tkerber
    a: int | 3 floats
95 1 tkerber
        The a-vector in scaled coordinates of the cell to cut out. If
96 1 tkerber
        integer, the a-vector will be the scaled vector from *origo* to the
97 1 tkerber
        atom with index *a*.
98 1 tkerber
    b: int | 3 floats
99 1 tkerber
        The b-vector in scaled coordinates of the cell to cut out. If
100 1 tkerber
        integer, the b-vector will be the scaled vector from *origo* to the
101 1 tkerber
        atom with index *b*.
102 1 tkerber
    c: int | 3 floats
103 1 tkerber
        The c-vector in scaled coordinates of the cell to cut out. If
104 1 tkerber
        integer, the c-vector will be the scaled vector from *origo* to the
105 1 tkerber
        atom with index *c*. Not used if *nlayers* is given.
106 1 tkerber
    origo: int | 3 floats
107 1 tkerber
        Position of origo of the new cell in scaled coordinates. If
108 1 tkerber
        integer, the position of the atom with index *origo* is used.
109 1 tkerber
    nlayers: int
110 1 tkerber
        If *nlayers* is not *None*, the returned cell will have
111 1 tkerber
        *nlayers* atomic layers in the c-direction. The direction of
112 1 tkerber
        the c-vector will be along cross(a, b) converted to real
113 1 tkerber
        space, i.e. normal to the plane spanned by a and b in
114 1 tkerber
        orthorombic systems.
115 1 tkerber
    extend: 1 or 3 floats
116 1 tkerber
        The *extend* argument scales the effective cell in which atoms
117 1 tkerber
        will be included. It must either be three floats or a single
118 1 tkerber
        float scaling all 3 directions.  By setting to a value just
119 1 tkerber
        above one, e.g. 1.05, it is possible to all the corner and
120 1 tkerber
        edge atoms in the returned cell.  This will of cause make the
121 1 tkerber
        returned cell non-repeatable, but is very usefull for
122 1 tkerber
        visualisation.
123 1 tkerber
    tolerance: float
124 1 tkerber
        Determines what is defined as a plane.  All atoms within
125 1 tkerber
        *tolerance* Angstroms from a given plane will be considered to
126 1 tkerber
        belong to that plane.
127 1 tkerber

128 1 tkerber
    Example:
129 1 tkerber

130 1 tkerber
    >>> import ase
131 1 tkerber
    >>> from ase.lattice.spacegroup import crystal
132 1 tkerber
    >>>
133 1 tkerber
    # Create an aluminium (111) slab with three layers
134 1 tkerber
    #
135 1 tkerber
    # First an unit cell of Al
136 1 tkerber
    >>> a = 4.05
137 1 tkerber
    >>> aluminium = crystal('Al', [(0,0,0)], spacegroup=225,
138 1 tkerber
    ...                     cellpar=[a, a, a, 90, 90, 90])
139 1 tkerber
    >>>
140 1 tkerber
    # Then cut out the slab
141 1 tkerber
    >>> al111 = cut(aluminium, (1,-1,0), (0,1,-1), nlayers=3)
142 1 tkerber
    >>>
143 1 tkerber
    # Visualisation of the skutterudite unit cell
144 1 tkerber
    #
145 1 tkerber
    # Again, create a skutterudite unit cell
146 1 tkerber
    >>> a = 9.04
147 1 tkerber
    >>> skutterudite = crystal(
148 1 tkerber
    ...     ('Co', 'Sb'),
149 1 tkerber
    ...     basis=[(0.25,0.25,0.25), (0.0, 0.335, 0.158)],
150 1 tkerber
    ...     spacegroup=204,
151 1 tkerber
    ...     cellpar=[a, a, a, 90, 90, 90])
152 1 tkerber
    >>>
153 1 tkerber
    # Then use *origo* to put 'Co' at the corners and *extend* to
154 1 tkerber
    # include all corner and edge atoms.
155 1 tkerber
    >>> s = cut(skutterudite, origo=(0.25, 0.25, 0.25), extend=1.01)
156 1 tkerber
    >>> ase.view(s)  # doctest: +SKIP
157 1 tkerber
    """
158 1 tkerber
    atoms = atoms.copy()
159 1 tkerber
    cell = atoms.cell
160 1 tkerber
161 1 tkerber
    if isinstance(origo, int):
162 1 tkerber
        origo = atoms.get_scaled_positions()[origo]
163 1 tkerber
    scaled = (atoms.get_scaled_positions() - origo)%1.0
164 1 tkerber
    scaled %= 1.0 # needed to ensure that all numbers are *less* than one
165 1 tkerber
    atoms.set_scaled_positions(scaled)
166 1 tkerber
167 1 tkerber
    if isinstance(a, int):
168 1 tkerber
        a = scaled[a] - origo
169 1 tkerber
    if isinstance(b, int):
170 1 tkerber
        b = scaled[b] - origo
171 1 tkerber
    if isinstance(c, int):
172 1 tkerber
        c = scaled[c] - origo
173 1 tkerber
174 1 tkerber
    a = np.array(a, dtype=float)
175 1 tkerber
    b = np.array(b, dtype=float)
176 1 tkerber
    origo = np.array(origo, dtype=float)
177 1 tkerber
178 1 tkerber
    if nlayers:
179 1 tkerber
        miller = np.cross(a, b) # surface normal
180 1 tkerber
        # The factor 36 = 2*2*3*3 is because the elements of a and b
181 1 tkerber
        # might be multiples of 1/2 or 1/3 because of lattice
182 1 tkerber
        # subtranslations
183 1 tkerber
        if np.all(36*miller - np.rint(36*miller)) < 1e-5:
184 1 tkerber
            miller = np.rint(36*miller)
185 1 tkerber
            miller /= gcd(miller)
186 1 tkerber
        tags, layers = get_layers(atoms, miller, tolerance)
187 1 tkerber
        while tags.max() < nlayers:
188 1 tkerber
            atoms = atoms.repeat(2)
189 1 tkerber
            tags, layers = get_layers(atoms, miller, tolerance)
190 1 tkerber
        # Convert surface normal in reciprocal space to direction in
191 1 tkerber
        # real space
192 1 tkerber
        metric = np.dot(cell, cell.T)
193 1 tkerber
        c = np.linalg.solve(metric.T, miller.T).T
194 1 tkerber
        c *= layers[nlayers]/np.sqrt(np.dot(c, miller))
195 1 tkerber
        if np.linalg.det(np.dot(np.array([a, b, c]), cell)) < 0:
196 1 tkerber
            c *= -1.0
197 1 tkerber
198 1 tkerber
    newcell = np.dot(np.array([a, b, c]), cell)
199 1 tkerber
200 1 tkerber
    # Create a new atoms object, repeated and translated such that
201 1 tkerber
    # it completely covers the new cell
202 1 tkerber
    scorners_newcell = np.array([[0., 0., 0.], [0., 0., 1.],
203 1 tkerber
                                 [0., 1., 0.], [0., 1., 1.],
204 1 tkerber
                                 [1., 0., 0.], [1., 0., 1.],
205 1 tkerber
                                 [1., 1., 0.], [1., 1., 1.]])
206 1 tkerber
    corners = np.dot(scorners_newcell, newcell*extend)
207 1 tkerber
    scorners = np.linalg.solve(cell.T, corners.T).T
208 1 tkerber
    rep = np.ceil(scorners.ptp(axis=0)).astype('int') + 1
209 1 tkerber
    trans = np.dot(np.floor(scorners.min(axis=0)), cell)
210 1 tkerber
    atoms = atoms.repeat(rep)
211 1 tkerber
    atoms.translate(trans)
212 1 tkerber
    atoms.set_cell(newcell)
213 1 tkerber
214 1 tkerber
    # Mask out atoms outside new cell
215 1 tkerber
    stol = tolerance  # scaled tolerance, XXX
216 1 tkerber
    maskcell = atoms.cell*extend
217 1 tkerber
    sp = np.linalg.solve(maskcell.T, (atoms.positions).T).T
218 1 tkerber
    mask = np.all(np.logical_and(-stol <= sp, sp < 1-stol), axis=1)
219 1 tkerber
    atoms = atoms[mask]
220 1 tkerber
221 1 tkerber
    return atoms
222 1 tkerber
223 1 tkerber
224 1 tkerber
225 1 tkerber
226 1 tkerber
def stack(atoms1, atoms2, axis=2, cell=None, fix=0.5,
227 1 tkerber
          maxstrain=0.5, distance=None):
228 1 tkerber
    """Return a new Atoms instance with *atoms2* added to atoms1
229 1 tkerber
    along the given axis. Periodicity in all directions is
230 1 tkerber
    ensured.
231 1 tkerber

232 1 tkerber
    The size of the final cell is determined by *cell*, except
233 1 tkerber
    that the length alongh *axis* will be the sum of
234 1 tkerber
    *atoms1.cell[axis]* and *atoms2.cell[axis]*. If *cell* is None,
235 1 tkerber
    it will be interpolated between *atoms1* and *atoms2*, where
236 1 tkerber
    *fix* determines their relative weight. Hence, if *fix* equals
237 1 tkerber
    zero, the final cell will be determined purely from *atoms1* and
238 1 tkerber
    if *fix* equals one, it will be determined purely from
239 1 tkerber
    *atoms2*.
240 1 tkerber

241 1 tkerber
    An ValueError exception will be raised if the far corner of
242 1 tkerber
    the unit cell of either *atoms1* or *atoms2* is displaced more
243 1 tkerber
    than *maxstrain*. Setting *maxstrain* to None, disable this
244 1 tkerber
    check.
245 1 tkerber

246 1 tkerber
    If *distance* is provided, the atomic positions in *atoms1* and
247 1 tkerber
    *atoms2* as well as the cell lengths along *axis* will be
248 1 tkerber
    adjusted such that the distance between the distance between
249 1 tkerber
    the closest atoms in *atoms1* and *atoms2* will equal *distance*.
250 1 tkerber
    This option uses scipy.optimize.fmin() and hence require scipy
251 1 tkerber
    to be installed.
252 1 tkerber

253 1 tkerber
    Example:
254 1 tkerber

255 1 tkerber
    >>> import ase
256 1 tkerber
    >>> from ase.lattice.spacegroup import crystal
257 1 tkerber
    >>>
258 1 tkerber
    # Create an Ag(110)-Si(110) interface with three atomic layers
259 1 tkerber
    # on each side.
260 1 tkerber
    >>> a_ag = 4.09
261 1 tkerber
    >>> ag = crystal(['Ag'], basis=[(0,0,0)], spacegroup=225,
262 1 tkerber
    ...              cellpar=[a_ag, a_ag, a_ag, 90., 90., 90.])
263 1 tkerber
    >>> ag110 = cut(ag, (0, 0, 3), (-1.5, 1.5, 0), nlayers=3)
264 1 tkerber
    >>>
265 1 tkerber
    >>> a_si = 5.43
266 1 tkerber
    >>> si = crystal(['Si'], basis=[(0,0,0)], spacegroup=227,
267 1 tkerber
    ...              cellpar=[a_si, a_si, a_si, 90., 90., 90.])
268 1 tkerber
    >>> si110 = cut(si, (0, 0, 2), (-1, 1, 0), nlayers=3)
269 1 tkerber
    >>>
270 1 tkerber
    >>> interface = stack(ag110, si110, maxstrain=1)
271 1 tkerber
    >>> ase.view(interface)  # doctest: +SKIP
272 1 tkerber
    >>>
273 1 tkerber
    # Once more, this time adjusted such that the distance between
274 1 tkerber
    # the closest Ag and Si atoms will be 2.3 Angstrom.
275 1 tkerber
    >>> interface2 = stack(ag110, si110,
276 1 tkerber
    ...                    maxstrain=1, distance=2.3)   # doctest:+ELLIPSIS
277 1 tkerber
    Optimization terminated successfully.
278 1 tkerber
        ...
279 1 tkerber
    >>> ase.view(interface2)  # doctest: +SKIP
280 1 tkerber
    """
281 1 tkerber
    atoms1 = atoms1.copy()
282 1 tkerber
    atoms2 = atoms2.copy()
283 1 tkerber
284 1 tkerber
    c1 = np.linalg.norm(atoms1.cell[axis])
285 1 tkerber
    c2 = np.linalg.norm(atoms2.cell[axis])
286 1 tkerber
    if cell is None:
287 1 tkerber
        cell1 = atoms1.cell.copy()
288 1 tkerber
        cell2 = atoms2.cell.copy()
289 1 tkerber
        cell1[axis] /= c1
290 1 tkerber
        cell2[axis] /= c2
291 1 tkerber
        cell = cell1 + fix*(cell2 - cell1)
292 1 tkerber
    cell[axis] /= np.linalg.norm(cell[axis])
293 1 tkerber
    cell1 = cell.copy()
294 1 tkerber
    cell2 = cell.copy()
295 1 tkerber
    cell1[axis] *= c1
296 1 tkerber
    cell2[axis] *= c2
297 1 tkerber
298 1 tkerber
    if (maxstrain and
299 1 tkerber
        (((cell1 - atoms1.cell).sum(axis=0)**2).sum() > maxstrain**2 or
300 1 tkerber
         ((cell2 - atoms2.cell).sum(axis=0)**2).sum() > maxstrain**2)):
301 1 tkerber
        raise ValueError('Incompatible cells.')
302 1 tkerber
303 1 tkerber
    sp1 = np.linalg.solve(atoms1.cell.T, atoms1.positions.T).T
304 1 tkerber
    sp2 = np.linalg.solve(atoms2.cell.T, atoms2.positions.T).T
305 1 tkerber
    atoms1.set_cell(cell1)
306 1 tkerber
    atoms2.set_cell(cell2)
307 1 tkerber
    atoms1.set_scaled_positions(sp1)
308 1 tkerber
    atoms2.set_scaled_positions(sp2)
309 1 tkerber
310 1 tkerber
    if distance is not None:
311 1 tkerber
        from scipy.optimize import fmin
312 1 tkerber
        def mindist(pos1, pos2):
313 1 tkerber
            n1 = len(pos1)
314 1 tkerber
            n2 = len(pos2)
315 1 tkerber
            idx1 = np.arange(n1).repeat(n2)
316 1 tkerber
            idx2 = np.tile(np.arange(n2), n1)
317 1 tkerber
            return np.sqrt(((pos1[idx1] - pos2[idx2])**2).sum(axis=1).min())
318 1 tkerber
        def func(x):
319 1 tkerber
            t1, t2, h1, h2 = x[0:3], x[3:6], x[6], x[7]
320 1 tkerber
            pos1 = atoms1.positions + t1
321 1 tkerber
            pos2 = atoms2.positions + t2
322 1 tkerber
            d1 = mindist(pos1, pos2 + (h1 + 1.0)*atoms1.cell[axis])
323 1 tkerber
            d2 = mindist(pos2, pos1 + (h2 + 1.0)*atoms2.cell[axis])
324 1 tkerber
            return (d1 - distance)**2 + (d2 - distance)**2
325 1 tkerber
        atoms1.center()
326 1 tkerber
        atoms2.center()
327 1 tkerber
        x0 = np.zeros((8,))
328 1 tkerber
        x = fmin(func, x0)
329 1 tkerber
        t1, t2, h1, h2 = x[0:3], x[3:6], x[6], x[7]
330 1 tkerber
        atoms1.translate(t1)
331 1 tkerber
        atoms2.translate(t2)
332 1 tkerber
        atoms1.cell[axis] *= 1.0 + h1
333 1 tkerber
        atoms2.cell[axis] *= 1.0 + h2
334 1 tkerber
335 1 tkerber
    atoms2.translate(atoms1.cell[axis])
336 1 tkerber
    atoms1.cell[axis] += atoms2.cell[axis]
337 1 tkerber
    atoms1.extend(atoms2)
338 1 tkerber
    return atoms1
339 1 tkerber
340 1 tkerber
341 1 tkerber
#-----------------------------------------------------------------
342 1 tkerber
# Self test
343 1 tkerber
if __name__ == '__main__':
344 1 tkerber
    import doctest
345 1 tkerber
    print 'doctest: ', doctest.testmod()