• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

materialsproject / pymatgen / 4075885785

pending completion
4075885785

push

github

Shyue Ping Ong
Merge branch 'master' of github.com:materialsproject/pymatgen

96 of 96 new or added lines in 27 files covered. (100.0%)

81013 of 102710 relevant lines covered (78.88%)

0.79 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

93.3
/pymatgen/core/surface.py
1
# Copyright (c) Pymatgen Development Team.
2
# Distributed under the terms of the MIT License.
3

4
"""
1✔
5
This module implements representations of slabs and surfaces, as well as
6
algorithms for generating them. If you use this module, please consider
7
citing the following work::
8

9
    R. Tran, Z. Xu, B. Radhakrishnan, D. Winston, W. Sun, K. A. Persson,
10
    S. P. Ong, "Surface Energies of Elemental Crystals", Scientific Data,
11
    2016, 3:160080, doi: 10.1038/sdata.2016.80.
12

13
as well as::
14

15
    Sun, W.; Ceder, G. Efficient creation and convergence of surface slabs,
16
    Surface Science, 2013, 617, 53-59, doi:10.1016/j.susc.2013.05.016.
17
"""
18

19
from __future__ import annotations
1✔
20

21
import copy
1✔
22
import itertools
1✔
23
import json
1✔
24
import logging
1✔
25
import math
1✔
26
import os
1✔
27
import warnings
1✔
28
from functools import reduce
1✔
29
from math import gcd
1✔
30

31
import numpy as np
1✔
32
from monty.fractions import lcm
1✔
33
from scipy.cluster.hierarchy import fcluster, linkage
1✔
34
from scipy.spatial.distance import squareform
1✔
35

36
from pymatgen.analysis.structure_matcher import StructureMatcher
1✔
37
from pymatgen.core.lattice import Lattice
1✔
38
from pymatgen.core.periodic_table import get_el_sp
1✔
39
from pymatgen.core.sites import PeriodicSite
1✔
40
from pymatgen.core.structure import Structure
1✔
41
from pymatgen.symmetry.analyzer import SpacegroupAnalyzer
1✔
42
from pymatgen.util.coord import in_coord_list
1✔
43

44
__author__ = "Richard Tran, Wenhao Sun, Zihan Xu, Shyue Ping Ong"
1✔
45

46

47
logger = logging.getLogger(__name__)
1✔
48

49

50
class Slab(Structure):
1✔
51
    """
52
    Subclass of Structure representing a Slab. Implements additional
53
    attributes pertaining to slabs, but the init method does not
54
    actually implement any algorithm that creates a slab. This is a
55
    DUMMY class who's init method only holds information about the
56
    slab. Also has additional methods that returns other information
57
    about a slab such as the surface area, normal, and atom adsorption.
58

59
    Note that all Slabs have the surface normal oriented perpendicular to the a
60
    and b lattice vectors. This means the lattice vectors a and b are in the
61
    surface plane and the c vector is out of the surface plane (though not
62
    necessarily perpendicular to the surface).
63

64
    .. attribute:: miller_index
65

66
        Miller index of plane parallel to surface.
67

68
    .. attribute:: scale_factor
69

70
        Final computed scale factor that brings the parent cell to the
71
        surface cell.
72

73
    .. attribute:: shift
74

75
        The shift value in Angstrom that indicates how much this
76
        slab has been shifted.
77
    """
78

79
    def __init__(
1✔
80
        self,
81
        lattice,
82
        species,
83
        coords,
84
        miller_index,
85
        oriented_unit_cell,
86
        shift,
87
        scale_factor,
88
        reorient_lattice=True,
89
        validate_proximity=False,
90
        to_unit_cell=False,
91
        reconstruction=None,
92
        coords_are_cartesian=False,
93
        site_properties=None,
94
        energy=None,
95
    ):
96
        """
97
        Makes a Slab structure, a structure object with additional information
98
        and methods pertaining to slabs.
99

100
        Args:
101
            lattice (Lattice/3x3 array): The lattice, either as a
102
                :class:`pymatgen.core.lattice.Lattice` or
103
                simply as any 2D array. Each row should correspond to a lattice
104
                vector. E.g., [[10,0,0], [20,10,0], [0,0,30]] specifies a
105
                lattice with lattice vectors [10,0,0], [20,10,0] and [0,0,30].
106
            species ([Species]): Sequence of species on each site. Can take in
107
                flexible input, including:
108

109
                i.  A sequence of element / species specified either as string
110
                    symbols, e.g. ["Li", "Fe2+", "P", ...] or atomic numbers,
111
                    e.g., (3, 56, ...) or actual Element or Species objects.
112

113
                ii. List of dict of elements/species and occupancies, e.g.,
114
                    [{"Fe" : 0.5, "Mn":0.5}, ...]. This allows the setup of
115
                    disordered structures.
116
            coords (Nx3 array): list of fractional/cartesian coordinates of
117
                each species.
118
            miller_index ([h, k, l]): Miller index of plane parallel to
119
                surface. Note that this is referenced to the input structure. If
120
                you need this to be based on the conventional cell,
121
                you should supply the conventional structure.
122
            oriented_unit_cell (Structure): The oriented_unit_cell from which
123
                this Slab is created (by scaling in the c-direction).
124
            shift (float): The shift in the c-direction applied to get the
125
                termination.
126
            scale_factor (np.ndarray): scale_factor Final computed scale factor
127
                that brings the parent cell to the surface cell.
128
            reorient_lattice (bool): reorients the lattice parameters such that
129
                the c direction is along the z axis.
130
            validate_proximity (bool): Whether to check if there are sites
131
                that are less than 0.01 Ang apart. Defaults to False.
132
            reconstruction (str): Type of reconstruction. Defaults to None if
133
                the slab is not reconstructed.
134
            coords_are_cartesian (bool): Set to True if you are providing
135
                coordinates in Cartesian coordinates. Defaults to False.
136
            site_properties (dict): Properties associated with the sites as a
137
                dict of sequences, e.g., {"magmom":[5,5,5,5]}. The sequences
138
                have to be the same length as the atomic species and
139
                fractional_coords. Defaults to None for no properties.
140
            energy (float): A value for the energy.
141
        """
142
        self.oriented_unit_cell = oriented_unit_cell
1✔
143
        self.miller_index = tuple(miller_index)
1✔
144
        self.shift = shift
1✔
145
        self.reconstruction = reconstruction
1✔
146
        self.scale_factor = np.array(scale_factor)
1✔
147
        self.energy = energy
1✔
148
        self.reorient_lattice = reorient_lattice
1✔
149
        if self.reorient_lattice:
1✔
150
            if coords_are_cartesian:
1✔
151
                coords = lattice.get_fractional_coords(coords)
1✔
152
                coords_are_cartesian = False
1✔
153
            lattice = Lattice.from_parameters(
1✔
154
                lattice.a,
155
                lattice.b,
156
                lattice.c,
157
                lattice.alpha,
158
                lattice.beta,
159
                lattice.gamma,
160
            )
161

162
        super().__init__(
1✔
163
            lattice,
164
            species,
165
            coords,
166
            validate_proximity=validate_proximity,
167
            to_unit_cell=to_unit_cell,
168
            coords_are_cartesian=coords_are_cartesian,
169
            site_properties=site_properties,
170
        )
171

172
    def get_orthogonal_c_slab(self):
1✔
173
        """
174
        This method returns a Slab where the normal (c lattice vector) is
175
        "forced" to be exactly orthogonal to the surface a and b lattice
176
        vectors. **Note that this breaks inherent symmetries in the slab.**
177
        It should be pointed out that orthogonality is not required to get good
178
        surface energies, but it can be useful in cases where the slabs are
179
        subsequently used for postprocessing of some kind, e.g. generating
180
        GBs or interfaces.
181
        """
182
        a, b, c = self.lattice.matrix
1✔
183
        new_c = np.cross(a, b)
1✔
184
        new_c /= np.linalg.norm(new_c)
1✔
185
        new_c = np.dot(c, new_c) * new_c
1✔
186
        new_latt = Lattice([a, b, new_c])
1✔
187
        return Slab(
1✔
188
            lattice=new_latt,
189
            species=self.species_and_occu,
190
            coords=self.cart_coords,
191
            miller_index=self.miller_index,
192
            oriented_unit_cell=self.oriented_unit_cell,
193
            shift=self.shift,
194
            scale_factor=self.scale_factor,
195
            coords_are_cartesian=True,
196
            energy=self.energy,
197
            reorient_lattice=self.reorient_lattice,
198
            site_properties=self.site_properties,
199
        )
200

201
    def get_tasker2_slabs(self, tol: float = 0.01, same_species_only=True):
1✔
202
        """
203
        Get a list of slabs that have been Tasker 2 corrected.
204

205
        Args:
206
            tol (float): Tolerance to determine if atoms are within same plane.
207
                This is a fractional tolerance, not an absolute one.
208
            same_species_only (bool): If True, only that are of the exact same
209
                species as the atom at the outermost surface are considered for
210
                moving. Otherwise, all atoms regardless of species that is
211
                within tol are considered for moving. Default is True (usually
212
                the desired behavior).
213

214
        Returns:
215
            ([Slab]) List of tasker 2 corrected slabs.
216
        """
217
        sites = list(self.sites)
1✔
218
        slabs = []
1✔
219

220
        sortedcsites = sorted(sites, key=lambda site: site.c)
1✔
221

222
        # Determine what fraction the slab is of the total cell size in the
223
        # c direction. Round to nearest rational number.
224
        nlayers_total = int(round(self.lattice.c / self.oriented_unit_cell.lattice.c))
1✔
225
        nlayers_slab = int(round((sortedcsites[-1].c - sortedcsites[0].c) * nlayers_total))
1✔
226
        slab_ratio = nlayers_slab / nlayers_total
1✔
227

228
        a = SpacegroupAnalyzer(self)
1✔
229
        symm_structure = a.get_symmetrized_structure()
1✔
230

231
        def equi_index(site):
1✔
232
            for i, equi_sites in enumerate(symm_structure.equivalent_sites):
1✔
233
                if site in equi_sites:
1✔
234
                    return i
1✔
235
            raise ValueError("Cannot determine equi index!")
×
236

237
        for surface_site, shift in [
1✔
238
            (sortedcsites[0], slab_ratio),
239
            (sortedcsites[-1], -slab_ratio),
240
        ]:
241
            tomove = []
1✔
242
            fixed = []
1✔
243
            for site in sites:
1✔
244
                if abs(site.c - surface_site.c) < tol and (
1✔
245
                    (not same_species_only) or site.species == surface_site.species
246
                ):
247
                    tomove.append(site)
1✔
248
                else:
249
                    fixed.append(site)
1✔
250

251
            # Sort and group the sites by the species and symmetry equivalence
252
            tomove = sorted(tomove, key=lambda s: equi_index(s))
1✔
253

254
            grouped = [list(sites) for k, sites in itertools.groupby(tomove, key=lambda s: equi_index(s))]
1✔
255

256
            if len(tomove) == 0 or any(len(g) % 2 != 0 for g in grouped):
1✔
257
                warnings.warn(
×
258
                    "Odd number of sites to divide! Try changing "
259
                    "the tolerance to ensure even division of "
260
                    "sites or create supercells in a or b directions "
261
                    "to allow for atoms to be moved!"
262
                )
263
                continue
×
264
            combinations = []
1✔
265
            for g in grouped:
1✔
266
                combinations.append(list(itertools.combinations(g, int(len(g) / 2))))
1✔
267

268
            for selection in itertools.product(*combinations):
1✔
269
                species = [site.species for site in fixed]
1✔
270
                fcoords = [site.frac_coords for site in fixed]
1✔
271

272
                for s in tomove:
1✔
273
                    species.append(s.species)
1✔
274
                    for group in selection:
1✔
275
                        if s in group:
1✔
276
                            fcoords.append(s.frac_coords)
1✔
277
                            break
1✔
278
                    else:
279
                        # Move unselected atom to the opposite surface.
280
                        fcoords.append(s.frac_coords + [0, 0, shift])
1✔
281

282
                # sort by species to put all similar species together.
283
                sp_fcoord = sorted(zip(species, fcoords), key=lambda x: x[0])
1✔
284
                species = [x[0] for x in sp_fcoord]
1✔
285
                fcoords = [x[1] for x in sp_fcoord]
1✔
286
                slab = Slab(
1✔
287
                    self.lattice,
288
                    species,
289
                    fcoords,
290
                    self.miller_index,
291
                    self.oriented_unit_cell,
292
                    self.shift,
293
                    self.scale_factor,
294
                    energy=self.energy,
295
                    reorient_lattice=self.reorient_lattice,
296
                )
297
                slabs.append(slab)
1✔
298
        s = StructureMatcher()
1✔
299
        unique = [ss[0] for ss in s.group_structures(slabs)]
1✔
300
        return unique
1✔
301

302
    def is_symmetric(self, symprec: float = 0.1):
1✔
303
        """
304
        Checks if surfaces are symmetric, i.e., contains inversion, mirror on (hkl) plane,
305
            or screw axis (rotation and translation) about [hkl].
306

307
        Args:
308
            symprec (float): Symmetry precision used for SpaceGroup analyzer.
309

310
        Returns:
311
            (bool) Whether surfaces are symmetric.
312
        """
313
        sg = SpacegroupAnalyzer(self, symprec=symprec)
1✔
314
        symmops = sg.get_point_group_operations()
1✔
315

316
        if (
1✔
317
            sg.is_laue()
318
            or any(op.translation_vector[2] != 0 for op in symmops)
319
            or any(np.alltrue(op.rotation_matrix[2] == np.array([0, 0, -1])) for op in symmops)
320
        ):
321
            # Check for inversion symmetry. Or if sites from surface (a) can be translated
322
            # to surface (b) along the [hkl]-axis, surfaces are symmetric. Or because the
323
            # two surfaces of our slabs are always parallel to the (hkl) plane,
324
            # any operation where there's an (hkl) mirror plane has surface symmetry
325
            return True
1✔
326
        return False
1✔
327

328
    def get_sorted_structure(self, key=None, reverse=False):
1✔
329
        """
330
        Get a sorted copy of the structure. The parameters have the same
331
        meaning as in list.sort. By default, sites are sorted by the
332
        electronegativity of the species. Note that Slab has to override this
333
        because of the different __init__ args.
334

335
        Args:
336
            key: Specifies a function of one argument that is used to extract
337
                a comparison key from each list element: key=str.lower. The
338
                default value is None (compare the elements directly).
339
            reverse (bool): If set to True, then the list elements are sorted
340
                as if each comparison were reversed.
341
        """
342
        sites = sorted(self, key=key, reverse=reverse)
1✔
343
        s = Structure.from_sites(sites)
1✔
344
        return Slab(
1✔
345
            s.lattice,
346
            s.species_and_occu,
347
            s.frac_coords,
348
            self.miller_index,
349
            self.oriented_unit_cell,
350
            self.shift,
351
            self.scale_factor,
352
            site_properties=s.site_properties,
353
            reorient_lattice=self.reorient_lattice,
354
        )
355

356
    def copy(self, site_properties=None, sanitize=False):
1✔
357
        """
358
        Convenience method to get a copy of the structure, with options to add
359
        site properties.
360

361
        Args:
362
            site_properties (dict): Properties to add or override. The
363
                properties are specified in the same way as the constructor,
364
                i.e., as a dict of the form {property: [values]}. The
365
                properties should be in the order of the *original* structure
366
                if you are performing sanitization.
367
            sanitize (bool): If True, this method will return a sanitized
368
                structure. Sanitization performs a few things: (i) The sites are
369
                sorted by electronegativity, (ii) a LLL lattice reduction is
370
                carried out to obtain a relatively orthogonalized cell,
371
                (iii) all fractional coords for sites are mapped into the
372
                unit cell.
373

374
        Returns:
375
            A copy of the Structure, with optionally new site_properties and
376
            optionally sanitized.
377
        """
378
        props = self.site_properties
1✔
379
        if site_properties:
1✔
380
            props.update(site_properties)
1✔
381
        return Slab(
1✔
382
            self.lattice,
383
            self.species_and_occu,
384
            self.frac_coords,
385
            self.miller_index,
386
            self.oriented_unit_cell,
387
            self.shift,
388
            self.scale_factor,
389
            site_properties=props,
390
            reorient_lattice=self.reorient_lattice,
391
        )
392

393
    @property
1✔
394
    def dipole(self):
1✔
395
        """
396
        Calculates the dipole of the Slab in the direction of the surface
397
        normal. Note that the Slab must be oxidation state-decorated for this
398
        to work properly. Otherwise, the Slab will always have a dipole of 0.
399
        """
400
        dipole = np.zeros(3)
1✔
401
        mid_pt = np.sum(self.cart_coords, axis=0) / len(self)
1✔
402
        normal = self.normal
1✔
403
        for site in self:
1✔
404
            charge = sum(getattr(sp, "oxi_state", 0) * amt for sp, amt in site.species.items())
1✔
405
            dipole += charge * np.dot(site.coords - mid_pt, normal) * normal
1✔
406
        return dipole
1✔
407

408
    def is_polar(self, tol_dipole_per_unit_area=1e-3) -> bool:
1✔
409
        """
410
        Checks whether the surface is polar by computing the dipole per unit
411
        area. Note that the Slab must be oxidation state-decorated for this
412
        to work properly. Otherwise, the Slab will always be non-polar.
413

414
        Args:
415
            tol_dipole_per_unit_area (float): A tolerance. If the dipole
416
                magnitude per unit area is less than this value, the Slab is
417
                considered non-polar. Defaults to 1e-3, which is usually
418
                pretty good. Normalized dipole per unit area is used as it is
419
                more reliable than using the total, which tends to be larger for
420
                slabs with larger surface areas.
421
        """
422
        dip_per_unit_area = self.dipole / self.surface_area
1✔
423
        return np.linalg.norm(dip_per_unit_area) > tol_dipole_per_unit_area
1✔
424

425
    @property
1✔
426
    def normal(self):
1✔
427
        """
428
        Calculates the surface normal vector of the slab
429
        """
430
        normal = np.cross(self.lattice.matrix[0], self.lattice.matrix[1])
1✔
431
        normal /= np.linalg.norm(normal)
1✔
432
        return normal
1✔
433

434
    @property
1✔
435
    def surface_area(self):
1✔
436
        """
437
        Calculates the surface area of the slab
438
        """
439
        m = self.lattice.matrix
1✔
440
        return np.linalg.norm(np.cross(m[0], m[1]))
1✔
441

442
    @property
1✔
443
    def center_of_mass(self):
1✔
444
        """
445
        Calculates the center of mass of the slab
446
        """
447
        weights = [s.species.weight for s in self]
1✔
448
        center_of_mass = np.average(self.frac_coords, weights=weights, axis=0)
1✔
449
        return center_of_mass
1✔
450

451
    def add_adsorbate_atom(self, indices, specie, distance):
1✔
452
        """
453
        Gets the structure of single atom adsorption.
454
        slab structure from the Slab class(in [0, 0, 1])
455

456
        Args:
457
            indices ([int]): Indices of sites on which to put the absorbate.
458
                Absorbed atom will be displaced relative to the center of
459
                these sites.
460
            specie (Species/Element/str): adsorbed atom species
461
            distance (float): between centers of the adsorbed atom and the
462
                given site in Angstroms.
463
        """
464
        # Let's do the work in Cartesian coords
465
        center = np.sum([self[i].coords for i in indices], axis=0) / len(indices)
1✔
466

467
        coords = center + self.normal * distance / np.linalg.norm(self.normal)
1✔
468

469
        self.append(specie, coords, coords_are_cartesian=True)
1✔
470

471
    def __str__(self):
1✔
472
        comp = self.composition
×
473
        outs = [
×
474
            f"Slab Summary ({comp.formula})",
475
            f"Reduced Formula: {comp.reduced_formula}",
476
            f"Miller index: {self.miller_index}",
477
            f"Shift: {self.shift:.4f}, Scale Factor: {self.scale_factor}",
478
        ]
479

480
        def to_s(x):
×
481
            return f"{x:0.6f}"
×
482

483
        outs.append("abc   : " + " ".join(to_s(i).rjust(10) for i in self.lattice.abc))
×
484
        outs.append("angles: " + " ".join(to_s(i).rjust(10) for i in self.lattice.angles))
×
485
        outs.append(f"Sites ({len(self)})")
×
486
        for i, site in enumerate(self):
×
487
            outs.append(
×
488
                " ".join(
489
                    [
490
                        str(i + 1),
491
                        site.species_string,
492
                        " ".join(to_s(j).rjust(12) for j in site.frac_coords),
493
                    ]
494
                )
495
            )
496
        return "\n".join(outs)
×
497

498
    def as_dict(self):
1✔
499
        """
500
        :return: MSONAble dict
501
        """
502
        d = super().as_dict()
1✔
503
        d["@module"] = type(self).__module__
1✔
504
        d["@class"] = type(self).__name__
1✔
505
        d["oriented_unit_cell"] = self.oriented_unit_cell.as_dict()
1✔
506
        d["miller_index"] = self.miller_index
1✔
507
        d["shift"] = self.shift
1✔
508
        d["scale_factor"] = self.scale_factor.tolist()
1✔
509
        d["reconstruction"] = self.reconstruction
1✔
510
        d["energy"] = self.energy
1✔
511
        return d
1✔
512

513
    @classmethod
1✔
514
    def from_dict(cls, d):
1✔
515
        """
516
        :param d: dict
517
        :return: Creates slab from dict.
518
        """
519
        lattice = Lattice.from_dict(d["lattice"])
1✔
520
        sites = [PeriodicSite.from_dict(sd, lattice) for sd in d["sites"]]
1✔
521
        s = Structure.from_sites(sites)
1✔
522

523
        return Slab(
1✔
524
            lattice=lattice,
525
            species=s.species_and_occu,
526
            coords=s.frac_coords,
527
            miller_index=d["miller_index"],
528
            oriented_unit_cell=Structure.from_dict(d["oriented_unit_cell"]),
529
            shift=d["shift"],
530
            scale_factor=d["scale_factor"],
531
            site_properties=s.site_properties,
532
            energy=d["energy"],
533
        )
534

535
    def get_surface_sites(self, tag=False):
1✔
536
        """
537
        Returns the surface sites and their indices in a dictionary. The
538
        oriented unit cell of the slab will determine the coordination number
539
        of a typical site. We use VoronoiNN to determine the
540
        coordination number of bulk sites and slab sites. Due to the
541
        pathological error resulting from some surface sites in the
542
        VoronoiNN, we assume any site that has this error is a surface
543
        site as well. This will work for elemental systems only for now. Useful
544
        for analysis involving broken bonds and for finding adsorption sites.
545

546
            Args:
547
                tag (bool): Option to adds site attribute "is_surfsite" (bool)
548
                    to all sites of slab. Defaults to False
549

550
            Returns:
551
                A dictionary grouping sites on top and bottom of the slab
552
                together.
553
                {"top": [sites with indices], "bottom": [sites with indices}
554

555
        TODO:
556
            Is there a way to determine site equivalence between sites in a slab
557
            and bulk system? This would allow us get the coordination number of
558
            a specific site for multi-elemental systems or systems with more
559
            than one unequivalent site. This will allow us to use this for
560
            compound systems.
561
        """
562
        from pymatgen.analysis.local_env import VoronoiNN
1✔
563

564
        # Get a dictionary of coordination numbers
565
        # for each distinct site in the structure
566
        a = SpacegroupAnalyzer(self.oriented_unit_cell)
1✔
567
        ucell = a.get_symmetrized_structure()
1✔
568
        cn_dict = {}
1✔
569
        v = VoronoiNN()
1✔
570
        unique_indices = [equ[0] for equ in ucell.equivalent_indices]
1✔
571

572
        for i in unique_indices:
1✔
573
            el = ucell[i].species_string
1✔
574
            if el not in cn_dict:
1✔
575
                cn_dict[el] = []
1✔
576
            # Since this will get the cn as a result of the weighted polyhedra, the
577
            # slightest difference in cn will indicate a different environment for a
578
            # species, eg. bond distance of each neighbor or neighbor species. The
579
            # decimal place to get some cn to be equal.
580
            cn = v.get_cn(ucell, i, use_weights=True)
1✔
581
            cn = float(f"{round(cn, 5):.5f}")
1✔
582
            if cn not in cn_dict[el]:
1✔
583
                cn_dict[el].append(cn)
1✔
584

585
        v = VoronoiNN()
1✔
586

587
        surf_sites_dict, properties = {"top": [], "bottom": []}, []
1✔
588
        for i, site in enumerate(self):
1✔
589
            # Determine if site is closer to the top or bottom of the slab
590
            top = site.frac_coords[2] > self.center_of_mass[2]
1✔
591

592
            try:
1✔
593
                # A site is a surface site, if its environment does
594
                # not fit the environment of other sites
595
                cn = float(f"{round(v.get_cn(self, i, use_weights=True), 5):.5f}")
1✔
596
                if cn < min(cn_dict[site.species_string]):
1✔
597
                    properties.append(True)
1✔
598
                    key = "top" if top else "bottom"
1✔
599
                    surf_sites_dict[key].append([site, i])
1✔
600
                else:
601
                    properties.append(False)
1✔
602
            except RuntimeError:
×
603
                # or if pathological error is returned, indicating a surface site
604
                properties.append(True)
×
605
                key = "top" if top else "bottom"
×
606
                surf_sites_dict[key].append([site, i])
×
607

608
        if tag:
1✔
609
            self.add_site_property("is_surf_site", properties)
×
610
        return surf_sites_dict
1✔
611

612
    def get_symmetric_site(self, point, cartesian=False):
1✔
613
        """
614
        This method uses symmetry operations to find equivalent sites on
615
            both sides of the slab. Works mainly for slabs with Laue
616
            symmetry. This is useful for retaining the non-polar and
617
            symmetric properties of a slab when creating adsorbed
618
            structures or symmetric reconstructions.
619

620
        Arg:
621
            point: Fractional coordinate.
622

623
        Returns:
624
            point: Fractional coordinate. A point equivalent to the
625
                parameter point, but on the other side of the slab
626
        """
627
        sg = SpacegroupAnalyzer(self)
1✔
628
        ops = sg.get_symmetry_operations(cartesian=cartesian)
1✔
629

630
        # Each operation on a point will return an equivalent point.
631
        # We want to find the point on the other side of the slab.
632
        for op in ops:
1✔
633
            slab = self.copy()
1✔
634
            site2 = op.operate(point)
1✔
635
            if f"{site2[2]:.6f}" == f"{point[2]:.6f}":
1✔
636
                continue
1✔
637

638
            # Add dummy site to check the overall structure is symmetric
639
            slab.append("O", point, coords_are_cartesian=cartesian)
1✔
640
            slab.append("O", site2, coords_are_cartesian=cartesian)
1✔
641
            sg = SpacegroupAnalyzer(slab)
1✔
642
            if sg.is_laue():
1✔
643
                break
1✔
644

645
            # If not symmetric, remove the two added
646
            # sites and try another symmetry operator
647
            slab.remove_sites([len(slab) - 1])
1✔
648
            slab.remove_sites([len(slab) - 1])
1✔
649

650
        return site2
1✔
651

652
    def symmetrically_add_atom(self, specie, point, coords_are_cartesian=False):
1✔
653
        """
654
        Class method for adding a site at a specified point in a slab.
655
            Will add the corresponding site on the other side of the
656
            slab to maintain equivalent surfaces.
657

658
        Arg:
659
            specie (str): The specie to add
660
            point (coords): The coordinate of the site in the slab to add.
661
            coords_are_cartesian (bool): Is the point in Cartesian coordinates
662

663
        Returns:
664
            (Slab): The modified slab
665
        """
666
        # For now just use the species of the
667
        # surface atom as the element to add
668

669
        # Get the index of the corresponding site at the bottom
670
        point2 = self.get_symmetric_site(point, cartesian=coords_are_cartesian)
1✔
671

672
        self.append(specie, point, coords_are_cartesian=coords_are_cartesian)
1✔
673
        self.append(specie, point2, coords_are_cartesian=coords_are_cartesian)
1✔
674

675
    def symmetrically_remove_atoms(self, indices):
1✔
676
        """
677
        Class method for removing sites corresponding to a list of indices.
678
            Will remove the corresponding site on the other side of the
679
            slab to maintain equivalent surfaces.
680

681
        Arg:
682
            indices ([indices]): The indices of the sites
683
                in the slab to remove.
684
        """
685
        slabcopy = SpacegroupAnalyzer(self.copy()).get_symmetrized_structure()
1✔
686
        points = [slabcopy[i].frac_coords for i in indices]
1✔
687
        removal_list = []
1✔
688

689
        for pt in points:
1✔
690
            # Get the index of the original site on top
691
            cart_point = slabcopy.lattice.get_cartesian_coords(pt)
1✔
692
            dist = [site.distance_from_point(cart_point) for site in slabcopy]
1✔
693
            site1 = dist.index(min(dist))
1✔
694

695
            # Get the index of the corresponding site at the bottom
696
            for i, eq_sites in enumerate(slabcopy.equivalent_sites):
1✔
697
                if slabcopy[site1] in eq_sites:
1✔
698
                    eq_indices = slabcopy.equivalent_indices[i]
1✔
699
                    break
1✔
700
            i1 = eq_indices[eq_sites.index(slabcopy[site1])]
1✔
701

702
            for i2 in eq_indices:
1✔
703
                if i2 == i1:
1✔
704
                    continue
×
705
                if slabcopy[i2].frac_coords[2] == slabcopy[i1].frac_coords[2]:
1✔
706
                    continue
×
707
                # Test site remove to see if it results in symmetric slab
708
                s = self.copy()
1✔
709
                s.remove_sites([i1, i2])
1✔
710
                if s.is_symmetric():
1✔
711
                    removal_list.extend([i1, i2])
1✔
712
                    break
1✔
713

714
        # If expected, 2 atoms are removed per index
715
        if len(removal_list) == 2 * len(indices):
1✔
716
            self.remove_sites(removal_list)
1✔
717
        else:
718
            warnings.warn("Equivalent sites could not be found for removal for all indices. Surface unchanged.")
×
719

720

721
class SlabGenerator:
1✔
722
    """
723
    This class generates different slabs using shift values determined by where
724
    a unique termination can be found along with other criteria such as where a
725
    termination doesn't break a polyhedral bond. The shift value then indicates
726
    where the slab layer will begin and terminate in the slab-vacuum system.
727

728
    .. attribute:: oriented_unit_cell
729

730
        A unit cell of the parent structure with the miller
731
        index of plane parallel to surface
732

733
    .. attribute:: parent
734

735
        Parent structure from which Slab was derived.
736

737
    .. attribute:: lll_reduce
738

739
        Whether or not the slabs will be orthogonalized
740

741
    .. attribute:: center_slab
742

743
        Whether or not the slabs will be centered between
744
        the vacuum layer
745

746
    .. attribute:: slab_scale_factor
747

748
        Final computed scale factor that brings the parent cell to the
749
        surface cell.
750

751
    .. attribute:: miller_index
752

753
        Miller index of plane parallel to surface.
754

755
    .. attribute:: min_slab_size
756

757
        Minimum size in angstroms of layers containing atoms
758

759
    .. attribute:: min_vac_size
760

761
        Minimize size in angstroms of layers containing vacuum
762

763
    """
764

765
    def __init__(
1✔
766
        self,
767
        initial_structure,
768
        miller_index,
769
        min_slab_size,
770
        min_vacuum_size,
771
        lll_reduce=False,
772
        center_slab=False,
773
        in_unit_planes=False,
774
        primitive=True,
775
        max_normal_search=None,
776
        reorient_lattice=True,
777
    ):
778
        """
779
        Calculates the slab scale factor and uses it to generate a unit cell
780
        of the initial structure that has been oriented by its miller index.
781
        Also stores the initial information needed later on to generate a slab.
782

783
        Args:
784
            initial_structure (Structure): Initial input structure. Note that to
785
                ensure that the miller indices correspond to usual
786
                crystallographic definitions, you should supply a conventional
787
                unit cell structure.
788
            miller_index ([h, k, l]): Miller index of plane parallel to
789
                surface. Note that this is referenced to the input structure. If
790
                you need this to be based on the conventional cell,
791
                you should supply the conventional structure.
792
            min_slab_size (float): In Angstroms or number of hkl planes
793
            min_vacuum_size (float): In Angstroms or number of hkl planes
794
            lll_reduce (bool): Whether to perform an LLL reduction on the
795
                eventual structure.
796
            center_slab (bool): Whether to center the slab in the cell with
797
                equal vacuum spacing from the top and bottom.
798
            in_unit_planes (bool): Whether to set min_slab_size and min_vac_size
799
                in units of hkl planes (True) or Angstrom (False/default).
800
                Setting in units of planes is useful for ensuring some slabs
801
                have a certain nlayer of atoms. e.g. for Cs (100), a 10 Ang
802
                slab will result in a slab with only 2 layer of atoms, whereas
803
                Fe (100) will have more layer of atoms. By using units of hkl
804
                planes instead, we ensure both slabs
805
                have the same number of atoms. The slab thickness will be in
806
                min_slab_size/math.ceil(self._proj_height/dhkl)
807
                multiples of oriented unit cells.
808
            primitive (bool): Whether to reduce any generated slabs to a
809
                primitive cell (this does **not** mean the slab is generated
810
                from a primitive cell, it simply means that after slab
811
                generation, we attempt to find shorter lattice vectors,
812
                which lead to less surface area and smaller cells).
813
            max_normal_search (int): If set to a positive integer, the code will
814
                conduct a search for a normal lattice vector that is as
815
                perpendicular to the surface as possible by considering
816
                multiples linear combinations of lattice vectors up to
817
                max_normal_search. This has no bearing on surface energies,
818
                but may be useful as a preliminary step to generating slabs
819
                for absorption and other sizes. It is typical that this will
820
                not be the smallest possible cell for simulation. Normality
821
                is not guaranteed, but the oriented cell will have the c
822
                vector as normal as possible (within the search range) to the
823
                surface. A value of up to the max absolute Miller index is
824
                usually sufficient.
825
            reorient_lattice (bool): reorients the lattice parameters such that
826
                the c direction is the third vector of the lattice matrix
827
        """
828
        # pylint: disable=E1130
829
        # Add Wyckoff symbols of the bulk, will help with
830
        # identfying types of sites in the slab system
831
        sg = SpacegroupAnalyzer(initial_structure)
1✔
832
        initial_structure.add_site_property("bulk_wyckoff", sg.get_symmetry_dataset()["wyckoffs"])
1✔
833
        initial_structure.add_site_property("bulk_equivalent", sg.get_symmetry_dataset()["equivalent_atoms"].tolist())
1✔
834
        latt = initial_structure.lattice
1✔
835
        miller_index = _reduce_vector(miller_index)
1✔
836
        # Calculate the surface normal using the reciprocal lattice vector.
837
        recp = latt.reciprocal_lattice_crystallographic
1✔
838
        normal = recp.get_cartesian_coords(miller_index)
1✔
839
        normal /= np.linalg.norm(normal)
1✔
840

841
        slab_scale_factor = []
1✔
842
        non_orth_ind = []
1✔
843
        eye = np.eye(3, dtype=int)
1✔
844
        for ii, jj in enumerate(miller_index):
1✔
845
            if jj == 0:
1✔
846
                # Lattice vector is perpendicular to surface normal, i.e.,
847
                # in plane of surface. We will simply choose this lattice
848
                # vector as one of the basis vectors.
849
                slab_scale_factor.append(eye[ii])
1✔
850
            else:
851
                # Calculate projection of lattice vector onto surface normal.
852
                d = abs(np.dot(normal, latt.matrix[ii])) / latt.abc[ii]
1✔
853
                non_orth_ind.append((ii, d))
1✔
854

855
        # We want the vector that has maximum magnitude in the
856
        # direction of the surface normal as the c-direction.
857
        # Results in a more "orthogonal" unit cell.
858
        c_index, dist = max(non_orth_ind, key=lambda t: t[1])
1✔
859

860
        if len(non_orth_ind) > 1:
1✔
861
            lcm_miller = lcm(*(miller_index[i] for i, d in non_orth_ind))
1✔
862
            for (ii, _di), (jj, _dj) in itertools.combinations(non_orth_ind, 2):
1✔
863
                l = [0, 0, 0]
1✔
864
                l[ii] = -int(round(lcm_miller / miller_index[ii]))
1✔
865
                l[jj] = int(round(lcm_miller / miller_index[jj]))
1✔
866
                slab_scale_factor.append(l)
1✔
867
                if len(slab_scale_factor) == 2:
1✔
868
                    break
1✔
869

870
        if max_normal_search is None:
1✔
871
            slab_scale_factor.append(eye[c_index])
1✔
872
        else:
873
            index_range = sorted(
1✔
874
                reversed(range(-max_normal_search, max_normal_search + 1)),
875
                key=lambda x: abs(x),
876
            )
877
            candidates = []
1✔
878
            for uvw in itertools.product(index_range, index_range, index_range):
1✔
879
                if (not any(uvw)) or abs(np.linalg.det(slab_scale_factor + [uvw])) < 1e-8:
1✔
880
                    continue
1✔
881
                vec = latt.get_cartesian_coords(uvw)
1✔
882
                l = np.linalg.norm(vec)
1✔
883
                cosine = abs(np.dot(vec, normal) / l)
1✔
884
                candidates.append((uvw, cosine, l))
1✔
885
                if abs(abs(cosine) - 1) < 1e-8:
1✔
886
                    # If cosine of 1 is found, no need to search further.
887
                    break
1✔
888
            # We want the indices with the maximum absolute cosine,
889
            # but smallest possible length.
890
            uvw, cosine, l = max(candidates, key=lambda x: (x[1], -x[2]))
1✔
891
            slab_scale_factor.append(uvw)
1✔
892

893
        slab_scale_factor = np.array(slab_scale_factor)
1✔
894

895
        # Let's make sure we have a left-handed crystallographic system
896
        if np.linalg.det(slab_scale_factor) < 0:
1✔
897
            slab_scale_factor *= -1
1✔
898

899
        # Make sure the slab_scale_factor is reduced to avoid
900
        # unnecessarily large slabs
901

902
        reduced_scale_factor = [_reduce_vector(v) for v in slab_scale_factor]
1✔
903
        slab_scale_factor = np.array(reduced_scale_factor)
1✔
904

905
        single = initial_structure.copy()
1✔
906
        single.make_supercell(slab_scale_factor)
1✔
907

908
        # When getting the OUC, lets return the most reduced
909
        # structure as possible to reduce calculations
910
        self.oriented_unit_cell = Structure.from_sites(single, to_unit_cell=True)
1✔
911
        self.max_normal_search = max_normal_search
1✔
912
        self.parent = initial_structure
1✔
913
        self.lll_reduce = lll_reduce
1✔
914
        self.center_slab = center_slab
1✔
915
        self.slab_scale_factor = slab_scale_factor
1✔
916
        self.miller_index = miller_index
1✔
917
        self.min_vac_size = min_vacuum_size
1✔
918
        self.min_slab_size = min_slab_size
1✔
919
        self.in_unit_planes = in_unit_planes
1✔
920
        self.primitive = primitive
1✔
921
        self._normal = normal
1✔
922
        a, b, c = self.oriented_unit_cell.lattice.matrix
1✔
923
        self._proj_height = abs(np.dot(normal, c))
1✔
924
        self.reorient_lattice = reorient_lattice
1✔
925

926
    def get_slab(self, shift=0, tol: float = 0.1, energy=None):
1✔
927
        """
928
        This method takes in shift value for the c lattice direction and
929
        generates a slab based on the given shift. You should rarely use this
930
        method. Instead, it is used by other generation algorithms to obtain
931
        all slabs.
932

933
        Arg:
934
            shift (float): A shift value in Angstrom that determines how much a
935
                slab should be shifted.
936
            tol (float): Tolerance to determine primitive cell.
937
            energy (float): An energy to assign to the slab.
938

939
        Returns:
940
            (Slab) A Slab object with a particular shifted oriented unit cell.
941
        """
942
        h = self._proj_height
1✔
943
        p = round(h / self.parent.lattice.d_hkl(self.miller_index), 8)
1✔
944
        if self.in_unit_planes:
1✔
945
            nlayers_slab = int(math.ceil(self.min_slab_size / p))
1✔
946
            nlayers_vac = int(math.ceil(self.min_vac_size / p))
1✔
947
        else:
948
            nlayers_slab = int(math.ceil(self.min_slab_size / h))
1✔
949
            nlayers_vac = int(math.ceil(self.min_vac_size / h))
1✔
950
        nlayers = nlayers_slab + nlayers_vac
1✔
951

952
        species = self.oriented_unit_cell.species_and_occu
1✔
953
        props = self.oriented_unit_cell.site_properties
1✔
954
        props = {k: v * nlayers_slab for k, v in props.items()}
1✔
955
        frac_coords = self.oriented_unit_cell.frac_coords
1✔
956
        frac_coords = np.array(frac_coords) + np.array([0, 0, -shift])[None, :]
1✔
957
        frac_coords -= np.floor(frac_coords)
1✔
958
        a, b, c = self.oriented_unit_cell.lattice.matrix
1✔
959
        new_lattice = [a, b, nlayers * c]
1✔
960
        frac_coords[:, 2] = frac_coords[:, 2] / nlayers
1✔
961
        all_coords = []
1✔
962
        for i in range(nlayers_slab):
1✔
963
            fcoords = frac_coords.copy()
1✔
964
            fcoords[:, 2] += i / nlayers
1✔
965
            all_coords.extend(fcoords)
1✔
966

967
        slab = Structure(new_lattice, species * nlayers_slab, all_coords, site_properties=props)
1✔
968

969
        scale_factor = self.slab_scale_factor
1✔
970
        # Whether or not to orthogonalize the structure
971
        if self.lll_reduce:
1✔
972
            lll_slab = slab.copy(sanitize=True)
×
973
            mapping = lll_slab.lattice.find_mapping(slab.lattice)
×
974
            scale_factor = np.dot(mapping[2], scale_factor)
×
975
            slab = lll_slab
×
976

977
        # Whether or not to center the slab layer around the vacuum
978
        if self.center_slab:
1✔
979
            avg_c = np.average([c[2] for c in slab.frac_coords])
1✔
980
            slab.translate_sites(list(range(len(slab))), [0, 0, 0.5 - avg_c])
1✔
981

982
        if self.primitive:
1✔
983
            prim = slab.get_primitive_structure(tolerance=tol)
1✔
984
            if energy is not None:
1✔
985
                energy = prim.volume / slab.volume * energy
1✔
986
            slab = prim
1✔
987

988
        # Reorient the lattice to get the correct reduced cell
989
        ouc = self.oriented_unit_cell.copy()
1✔
990
        if self.primitive:
1✔
991
            # find a reduced ouc
992
            slab_l = slab.lattice
1✔
993
            ouc = ouc.get_primitive_structure(
1✔
994
                constrain_latt={
995
                    "a": slab_l.a,
996
                    "b": slab_l.b,
997
                    "alpha": slab_l.alpha,
998
                    "beta": slab_l.beta,
999
                    "gamma": slab_l.gamma,
1000
                }
1001
            )
1002
            # Check this is the correct oriented unit cell
1003
            ouc = self.oriented_unit_cell if slab_l.a != ouc.lattice.a or slab_l.b != ouc.lattice.b else ouc
1✔
1004

1005
        return Slab(
1✔
1006
            slab.lattice,
1007
            slab.species_and_occu,
1008
            slab.frac_coords,
1009
            self.miller_index,
1010
            ouc,
1011
            shift,
1012
            scale_factor,
1013
            energy=energy,
1014
            site_properties=slab.site_properties,
1015
            reorient_lattice=self.reorient_lattice,
1016
        )
1017

1018
    def _calculate_possible_shifts(self, tol: float = 0.1):
1✔
1019
        frac_coords = self.oriented_unit_cell.frac_coords
1✔
1020
        n = len(frac_coords)
1✔
1021

1022
        if n == 1:
1✔
1023
            # Clustering does not work when there is only one data point.
1024
            shift = frac_coords[0][2] + 0.5
1✔
1025
            return [shift - math.floor(shift)]
1✔
1026

1027
        # We cluster the sites according to the c coordinates. But we need to
1028
        # take into account PBC. Let's compute a fractional c-coordinate
1029
        # distance matrix that accounts for PBC.
1030
        dist_matrix = np.zeros((n, n))
1✔
1031
        h = self._proj_height
1✔
1032
        # Projection of c lattice vector in
1033
        # direction of surface normal.
1034
        for i, j in itertools.combinations(list(range(n)), 2):
1✔
1035
            if i != j:
1✔
1036
                cdist = frac_coords[i][2] - frac_coords[j][2]
1✔
1037
                cdist = abs(cdist - round(cdist)) * h
1✔
1038
                dist_matrix[i, j] = cdist
1✔
1039
                dist_matrix[j, i] = cdist
1✔
1040

1041
        condensed_m = squareform(dist_matrix)
1✔
1042
        z = linkage(condensed_m)
1✔
1043
        clusters = fcluster(z, tol, criterion="distance")
1✔
1044

1045
        # Generate dict of cluster# to c val - doesn't matter what the c is.
1046
        c_loc = {c: frac_coords[i][2] for i, c in enumerate(clusters)}
1✔
1047

1048
        # Put all c into the unit cell.
1049
        possible_c = [c - math.floor(c) for c in sorted(c_loc.values())]
1✔
1050

1051
        # Calculate the shifts
1052
        nshifts = len(possible_c)
1✔
1053
        shifts = []
1✔
1054
        for i in range(nshifts):
1✔
1055
            if i == nshifts - 1:
1✔
1056
                # There is an additional shift between the first and last c
1057
                # coordinate. But this needs special handling because of PBC.
1058
                shift = (possible_c[0] + 1 + possible_c[i]) * 0.5
1✔
1059
                if shift > 1:
1✔
1060
                    shift -= 1
1✔
1061
            else:
1062
                shift = (possible_c[i] + possible_c[i + 1]) * 0.5
1✔
1063
            shifts.append(shift - math.floor(shift))
1✔
1064
        shifts = sorted(shifts)
1✔
1065
        return shifts
1✔
1066

1067
    def _get_c_ranges(self, bonds):
1✔
1068
        c_ranges = []
1✔
1069
        bonds = {(get_el_sp(s1), get_el_sp(s2)): dist for (s1, s2), dist in bonds.items()}
1✔
1070
        for (sp1, sp2), bond_dist in bonds.items():
1✔
1071
            for site in self.oriented_unit_cell:
1✔
1072
                if sp1 in site.species:
1✔
1073
                    for nn in self.oriented_unit_cell.get_neighbors(site, bond_dist):
1✔
1074
                        if sp2 in nn.species:
1✔
1075
                            c_range = tuple(sorted([site.frac_coords[2], nn.frac_coords[2]]))
1✔
1076
                            if c_range[1] > 1:
1✔
1077
                                # Takes care of PBC when c coordinate of site
1078
                                # goes beyond the upper boundary of the cell
1079
                                c_ranges.append((c_range[0], 1))
1✔
1080
                                c_ranges.append((0, c_range[1] - 1))
1✔
1081
                            elif c_range[0] < 0:
1✔
1082
                                # Takes care of PBC when c coordinate of site
1083
                                # is below the lower boundary of the unit cell
1084
                                c_ranges.append((0, c_range[1]))
1✔
1085
                                c_ranges.append((c_range[0] + 1, 1))
1✔
1086
                            elif c_range[0] != c_range[1]:
1✔
1087
                                c_ranges.append((c_range[0], c_range[1]))
1✔
1088
        return c_ranges
1✔
1089

1090
    def get_slabs(
1✔
1091
        self,
1092
        bonds=None,
1093
        ftol=0.1,
1094
        tol=0.1,
1095
        max_broken_bonds=0,
1096
        symmetrize=False,
1097
        repair=False,
1098
    ):
1099
        """
1100
        This method returns a list of slabs that are generated using the list of
1101
        shift values from the method, _calculate_possible_shifts(). Before the
1102
        shifts are used to create the slabs however, if the user decides to take
1103
        into account whether or not a termination will break any polyhedral
1104
        structure (bonds is not None), this method will filter out any shift
1105
        values that do so.
1106

1107
        Args:
1108
            bonds ({(specie1, specie2): max_bond_dist}: bonds are
1109
                specified as a dict of tuples: float of specie1, specie2
1110
                and the max bonding distance. For example, PO4 groups may be
1111
                defined as {("P", "O"): 3}.
1112
            tol (float): General tolerance parameter for getting primitive
1113
                cells and matching structures
1114
            ftol (float): Threshold parameter in fcluster in order to check
1115
                if two atoms are lying on the same plane. Default thresh set
1116
                to 0.1 Angstrom in the direction of the surface normal.
1117
            max_broken_bonds (int): Maximum number of allowable broken bonds
1118
                for the slab. Use this to limit # of slabs (some structures
1119
                may have a lot of slabs). Defaults to zero, which means no
1120
                defined bonds must be broken.
1121
            symmetrize (bool): Whether or not to ensure the surfaces of the
1122
                slabs are equivalent.
1123
            repair (bool): Whether to repair terminations with broken bonds
1124
                or just omit them. Set to False as repairing terminations can
1125
                lead to many possible slabs as oppose to just omitting them.
1126

1127
        Returns:
1128
            ([Slab]) List of all possible terminations of a particular surface.
1129
            Slabs are sorted by the # of bonds broken.
1130
        """
1131
        c_ranges = [] if bonds is None else self._get_c_ranges(bonds)
1✔
1132

1133
        slabs = []
1✔
1134
        for shift in self._calculate_possible_shifts(tol=ftol):
1✔
1135
            bonds_broken = 0
1✔
1136
            for r in c_ranges:
1✔
1137
                if r[0] <= shift <= r[1]:
1✔
1138
                    bonds_broken += 1
1✔
1139
            slab = self.get_slab(shift, tol=tol, energy=bonds_broken)
1✔
1140
            if bonds_broken <= max_broken_bonds:
1✔
1141
                slabs.append(slab)
1✔
1142
            elif repair:
1✔
1143
                # If the number of broken bonds is exceeded,
1144
                # we repair the broken bonds on the slab
1145
                slabs.append(self.repair_broken_bonds(slab, bonds))
1✔
1146

1147
        # Further filters out any surfaces made that might be the same
1148
        m = StructureMatcher(ltol=tol, stol=tol, primitive_cell=False, scale=False)
1✔
1149

1150
        new_slabs = []
1✔
1151
        for g in m.group_structures(slabs):
1✔
1152
            # For each unique termination, symmetrize the
1153
            # surfaces by removing sites from the bottom.
1154
            if symmetrize:
1✔
1155
                slabs = self.nonstoichiometric_symmetrized_slab(g[0])
1✔
1156
                new_slabs.extend(slabs)
1✔
1157
            else:
1158
                new_slabs.append(g[0])
1✔
1159

1160
        match = StructureMatcher(ltol=tol, stol=tol, primitive_cell=False, scale=False)
1✔
1161
        new_slabs = [g[0] for g in match.group_structures(new_slabs)]
1✔
1162

1163
        return sorted(new_slabs, key=lambda s: s.energy)
1✔
1164

1165
    def repair_broken_bonds(self, slab, bonds):
1✔
1166
        """
1167
        This method will find undercoordinated atoms due to slab
1168
        cleaving specified by the bonds parameter and move them
1169
        to the other surface to make sure the bond is kept intact.
1170
        In a future release of surface.py, the ghost_sites will be
1171
        used to tell us how the repair bonds should look like.
1172

1173
        Arg:
1174
            slab (structure): A structure object representing a slab.
1175
            bonds ({(specie1, specie2): max_bond_dist}: bonds are
1176
                specified as a dict of tuples: float of specie1, specie2
1177
                and the max bonding distance. For example, PO4 groups may be
1178
                defined as {("P", "O"): 3}.
1179

1180
        Returns:
1181
            (Slab) A Slab object with a particular shifted oriented unit cell.
1182
        """
1183
        for pair in bonds:
1✔
1184
            blength = bonds[pair]
1✔
1185

1186
            # First lets determine which element should be the
1187
            # reference (center element) to determine broken bonds.
1188
            # e.g. P for a PO4 bond. Find integer coordination
1189
            # numbers of the pair of elements wrt to each other
1190
            cn_dict = {}
1✔
1191
            for i, el in enumerate(pair):
1✔
1192
                cnlist = []
1✔
1193
                for site in self.oriented_unit_cell:
1✔
1194
                    poly_coord = 0
1✔
1195
                    if site.species_string == el:
1✔
1196
                        for nn in self.oriented_unit_cell.get_neighbors(site, blength):
1✔
1197
                            if nn[0].species_string == pair[i - 1]:
1✔
1198
                                poly_coord += 1
1✔
1199
                    cnlist.append(poly_coord)
1✔
1200
                cn_dict[el] = cnlist
1✔
1201

1202
            # We make the element with the higher coordination our reference
1203
            if max(cn_dict[pair[0]]) > max(cn_dict[pair[1]]):
1✔
1204
                element1, element2 = pair
1✔
1205
            else:
1206
                element2, element1 = pair
×
1207

1208
            for i, site in enumerate(slab):
1✔
1209
                # Determine the coordination of our reference
1210
                if site.species_string == element1:
1✔
1211
                    poly_coord = 0
1✔
1212
                    for neighbor in slab.get_neighbors(site, blength):
1✔
1213
                        poly_coord += 1 if neighbor.species_string == element2 else 0
1✔
1214

1215
                    # suppose we find an undercoordinated reference atom
1216
                    if poly_coord not in cn_dict[element1]:
1✔
1217
                        # We get the reference atom of the broken bonds
1218
                        # (undercoordinated), move it to the other surface
1219
                        slab = self.move_to_other_side(slab, [i])
1✔
1220

1221
                        # find its NNs with the corresponding
1222
                        # species it should be coordinated with
1223
                        neighbors = slab.get_neighbors(slab[i], blength, include_index=True)
1✔
1224
                        tomove = [nn[2] for nn in neighbors if nn[0].species_string == element2]
1✔
1225
                        tomove.append(i)
1✔
1226
                        # and then move those NNs along with the central
1227
                        # atom back to the other side of the slab again
1228
                        slab = self.move_to_other_side(slab, tomove)
1✔
1229

1230
        return slab
1✔
1231

1232
    def move_to_other_side(self, init_slab, index_of_sites):
1✔
1233
        """
1234
        This method will Move a set of sites to the
1235
        other side of the slab (opposite surface).
1236

1237
        Arg:
1238
            init_slab (structure): A structure object representing a slab.
1239
            index_of_sites (list of ints): The list of indices representing
1240
                the sites we want to move to the other side.
1241

1242
        Returns:
1243
            (Slab) A Slab object with a particular shifted oriented unit cell.
1244
        """
1245
        slab = init_slab.copy()
1✔
1246

1247
        # Determine what fraction the slab is of the total cell size
1248
        # in the c direction. Round to nearest rational number.
1249
        h = self._proj_height
1✔
1250
        p = h / self.parent.lattice.d_hkl(self.miller_index)
1✔
1251
        if self.in_unit_planes:
1✔
1252
            nlayers_slab = int(math.ceil(self.min_slab_size / p))
×
1253
            nlayers_vac = int(math.ceil(self.min_vac_size / p))
×
1254
        else:
1255
            nlayers_slab = int(math.ceil(self.min_slab_size / h))
1✔
1256
            nlayers_vac = int(math.ceil(self.min_vac_size / h))
1✔
1257
        nlayers = nlayers_slab + nlayers_vac
1✔
1258
        slab_ratio = nlayers_slab / nlayers
1✔
1259

1260
        # Sort the index of sites based on which side they are on
1261
        top_site_index = [i for i in index_of_sites if slab[i].frac_coords[2] > slab.center_of_mass[2]]
1✔
1262
        bottom_site_index = [i for i in index_of_sites if slab[i].frac_coords[2] < slab.center_of_mass[2]]
1✔
1263

1264
        # Translate sites to the opposite surfaces
1265
        slab.translate_sites(top_site_index, [0, 0, slab_ratio])
1✔
1266
        slab.translate_sites(bottom_site_index, [0, 0, -slab_ratio])
1✔
1267

1268
        return Slab(
1✔
1269
            init_slab.lattice,
1270
            slab.species,
1271
            slab.frac_coords,
1272
            init_slab.miller_index,
1273
            init_slab.oriented_unit_cell,
1274
            init_slab.shift,
1275
            init_slab.scale_factor,
1276
            energy=init_slab.energy,
1277
        )
1278

1279
    def nonstoichiometric_symmetrized_slab(self, init_slab):
1✔
1280
        """
1281
        This method checks whether or not the two surfaces of the slab are
1282
        equivalent. If the point group of the slab has an inversion symmetry (
1283
        ie. belong to one of the Laue groups), then it is assumed that the
1284
        surfaces should be equivalent. Otherwise, sites at the bottom of the
1285
        slab will be removed until the slab is symmetric. Note the removal of sites
1286
        can destroy the stoichiometry of the slab. For non-elemental
1287
        structures, the chemical potential will be needed to calculate surface energy.
1288

1289
        Arg:
1290
            init_slab (Structure): A single slab structure
1291

1292
        Returns:
1293
            Slab (structure): A symmetrized Slab object.
1294
        """
1295
        if init_slab.is_symmetric():
1✔
1296
            return [init_slab]
1✔
1297

1298
        nonstoich_slabs = []
1✔
1299
        # Build an equivalent surface slab for each of the different surfaces
1300
        for top in [True, False]:
1✔
1301
            asym = True
1✔
1302
            slab = init_slab.copy()
1✔
1303
            slab.energy = init_slab.energy
1✔
1304

1305
            while asym:
1✔
1306
                # Keep removing sites from the bottom one by one until both
1307
                # surfaces are symmetric or the number of sites removed has
1308
                # exceeded 10 percent of the original slab
1309

1310
                c_dir = [site[2] for i, site in enumerate(slab.frac_coords)]
1✔
1311

1312
                if top:
1✔
1313
                    slab.remove_sites([c_dir.index(max(c_dir))])
1✔
1314
                else:
1315
                    slab.remove_sites([c_dir.index(min(c_dir))])
1✔
1316
                if len(slab) <= len(self.parent):
1✔
1317
                    break
1✔
1318

1319
                # Check if the altered surface is symmetric
1320
                if slab.is_symmetric():
1✔
1321
                    asym = False
1✔
1322
                    nonstoich_slabs.append(slab)
1✔
1323

1324
        if len(slab) <= len(self.parent):
1✔
1325
            warnings.warn("Too many sites removed, please use a larger slab size.")
1✔
1326

1327
        return nonstoich_slabs
1✔
1328

1329

1330
module_dir = os.path.dirname(os.path.abspath(__file__))
1✔
1331
with open(os.path.join(module_dir, "reconstructions_archive.json")) as data_file:
1✔
1332
    reconstructions_archive = json.load(data_file)
1✔
1333

1334

1335
class ReconstructionGenerator:
1✔
1336
    """
1337
    This class takes in a pre-defined dictionary specifying the parameters
1338
    need to build a reconstructed slab such as the SlabGenerator parameters,
1339
    transformation matrix, sites to remove/add and slab/vacuum size. It will
1340
    then use the formatted instructions provided by the dictionary to build
1341
    the desired reconstructed slab from the initial structure.
1342

1343
    .. attribute:: slabgen_params
1344

1345
        Parameters for the SlabGenerator
1346

1347
    .. trans_matrix::
1348

1349
        A 3x3 transformation matrix to generate the reconstructed
1350
            slab. Only the a and b lattice vectors are actually
1351
            changed while the c vector remains the same. This
1352
            matrix is what the Wood's notation is based on.
1353

1354
    .. reconstruction_json::
1355

1356
        The full json or dictionary containing the instructions
1357
            for building the reconstructed slab
1358

1359
    .. termination::
1360

1361
        The index of the termination of the slab
1362

1363
    TODO:
1364
    - Right now there is no way to specify what atom is being
1365
        added. In the future, use basis sets?
1366
    """
1367

1368
    def __init__(self, initial_structure, min_slab_size, min_vacuum_size, reconstruction_name):
1✔
1369
        """
1370
        Generates reconstructed slabs from a set of instructions
1371
            specified by a dictionary or json file.
1372

1373
        Args:
1374
            initial_structure (Structure): Initial input structure. Note
1375
                that to ensure that the miller indices correspond to usual
1376
                crystallographic definitions, you should supply a conventional
1377
                unit cell structure.
1378
            min_slab_size (float): In Angstroms
1379
            min_vacuum_size (float): In Angstroms
1380

1381
            reconstruction (str): Name of the dict containing the instructions
1382
                for building a reconstructed slab. The dictionary can contain
1383
                any item the creator deems relevant, however any instructions
1384
                archived in pymatgen for public use needs to contain the
1385
                following keys and items to ensure compatibility with the
1386
                ReconstructionGenerator:
1387

1388
                    "name" (str): A descriptive name for the type of
1389
                        reconstruction. Typically the name will have the type
1390
                        of structure the reconstruction is for, the Miller
1391
                        index, and Wood's notation along with anything to
1392
                        describe the reconstruction: e.g.:
1393
                        "fcc_110_missing_row_1x2"
1394
                    "description" (str): A longer description of your
1395
                        reconstruction. This is to help future contributors who
1396
                        want to add other types of reconstructions to the
1397
                        archive on pymatgen to check if the reconstruction
1398
                        already exists. Please read the descriptions carefully
1399
                        before adding a new type of reconstruction to ensure it
1400
                        is not in the archive yet.
1401
                    "reference" (str): Optional reference to where the
1402
                        reconstruction was taken from or first observed.
1403
                    "spacegroup" (dict): e.g. {"symbol": "Fm-3m", "number": 225}
1404
                        Indicates what kind of structure is this reconstruction.
1405
                    "miller_index" ([h,k,l]): Miller index of your reconstruction
1406
                    "Woods_notation" (str): For a reconstruction, the a and b
1407
                        lattice may change to accommodate the symmetry of the
1408
                        reconstruction. This notation indicates the change in
1409
                        the vectors relative to the primitive (p) or
1410
                        conventional (c) slab cell. E.g. p(2x1):
1411

1412
                        Wood, E. A. (1964). Vocabulary of surface
1413
                        crystallography. Journal of Applied Physics, 35(4),
1414
                        1306-1312.
1415

1416
                    "transformation_matrix" (numpy array): A 3x3 matrix to
1417
                        transform the slab. Only the a and b lattice vectors
1418
                        should change while the c vector remains the same.
1419
                    "SlabGenerator_parameters" (dict): A dictionary containing
1420
                        the parameters for the SlabGenerator class excluding the
1421
                        miller_index, min_slab_size and min_vac_size as the
1422
                        Miller index is already specified and the min_slab_size
1423
                        and min_vac_size can be changed regardless of what type
1424
                        of reconstruction is used. Having a consistent set of
1425
                        SlabGenerator parameters allows for the instructions to
1426
                        be reused to consistently build a reconstructed slab.
1427
                    "points_to_remove" (list of coords): A list of sites to
1428
                        remove where the first two indices are fraction (in a
1429
                        and b) and the third index is in units of 1/d (in c).
1430
                    "points_to_add" (list of frac_coords): A list of sites to add
1431
                        where the first two indices are fraction (in a an b) and
1432
                        the third index is in units of 1/d (in c).
1433

1434
                    "base_reconstruction" (dict): Option to base a reconstruction on
1435
                        an existing reconstruction model also exists to easily build
1436
                        the instructions without repeating previous work. E.g. the
1437
                        alpha reconstruction of halites is based on the octopolar
1438
                        reconstruction but with the topmost atom removed. The dictionary
1439
                        for the alpha reconstruction would therefore contain the item
1440
                        "reconstruction_base": "halite_111_octopolar_2x2", and
1441
                        additional sites for "points_to_remove" and "points_to_add"
1442
                        can be added to modify this reconstruction.
1443

1444
                    For "points_to_remove" and "points_to_add", the third index for
1445
                        the c vector is in units of 1/d where d is the spacing
1446
                        between atoms along hkl (the c vector) and is relative to
1447
                        the topmost site in the unreconstructed slab. e.g. a point
1448
                        of [0.5, 0.25, 1] corresponds to the 0.5 frac_coord of a,
1449
                        0.25 frac_coord of b and a distance of 1 atomic layer above
1450
                        the topmost site. [0.5, 0.25, -0.5] where the third index
1451
                        corresponds to a point half a atomic layer below the topmost
1452
                        site. [0.5, 0.25, 0] corresponds to a point in the same
1453
                        position along c as the topmost site. This is done because
1454
                        while the primitive units of a and b will remain constant,
1455
                        the user can vary the length of the c direction by changing
1456
                        the slab layer or the vacuum layer.
1457

1458
            NOTE: THE DICTIONARY SHOULD ONLY CONTAIN "points_to_remove" AND
1459
            "points_to_add" FOR THE TOP SURFACE. THE ReconstructionGenerator
1460
            WILL MODIFY THE BOTTOM SURFACE ACCORDINGLY TO RETURN A SLAB WITH
1461
            EQUIVALENT SURFACES.
1462
        """
1463
        if reconstruction_name not in reconstructions_archive:
1✔
1464
            raise KeyError(
×
1465
                f"The reconstruction_name entered ({reconstruction_name}) does not exist in the "
1466
                f"archive. Please select from one of the following reconstructions: {list(reconstructions_archive)} "
1467
                "or add the appropriate dictionary to the archive file "
1468
                "reconstructions_archive.json."
1469
            )
1470

1471
        # Get the instructions to build the reconstruction
1472
        # from the reconstruction_archive
1473
        recon_json = copy.deepcopy(reconstructions_archive[reconstruction_name])
1✔
1474
        new_points_to_add, new_points_to_remove = [], []
1✔
1475
        if "base_reconstruction" in recon_json:
1✔
1476
            if "points_to_add" in recon_json:
1✔
1477
                new_points_to_add = recon_json["points_to_add"]
1✔
1478
            if "points_to_remove" in recon_json:
1✔
1479
                new_points_to_remove = recon_json["points_to_remove"]
×
1480

1481
            # Build new instructions from a base reconstruction
1482
            recon_json = copy.deepcopy(reconstructions_archive[recon_json["base_reconstruction"]])
1✔
1483
            if "points_to_add" in recon_json:
1✔
1484
                del recon_json["points_to_add"]
1✔
1485
            if "points_to_remove" in recon_json:
1✔
1486
                del recon_json["points_to_remove"]
×
1487
            if new_points_to_add:
1✔
1488
                recon_json["points_to_add"] = new_points_to_add
1✔
1489
            if new_points_to_remove:
1✔
1490
                recon_json["points_to_remove"] = new_points_to_remove
×
1491

1492
        slabgen_params = copy.deepcopy(recon_json["SlabGenerator_parameters"])
1✔
1493
        slabgen_params["initial_structure"] = initial_structure.copy()
1✔
1494
        slabgen_params["miller_index"] = recon_json["miller_index"]
1✔
1495
        slabgen_params["min_slab_size"] = min_slab_size
1✔
1496
        slabgen_params["min_vacuum_size"] = min_vacuum_size
1✔
1497

1498
        self.slabgen_params = slabgen_params
1✔
1499
        self.trans_matrix = recon_json["transformation_matrix"]
1✔
1500
        self.reconstruction_json = recon_json
1✔
1501
        self.name = reconstruction_name
1✔
1502

1503
    def build_slabs(self):
1✔
1504
        """
1505
        Builds the reconstructed slab by:
1506
            (1) Obtaining the unreconstructed slab using the specified
1507
                parameters for the SlabGenerator.
1508
            (2) Applying the appropriate lattice transformation in the
1509
                a and b lattice vectors.
1510
            (3) Remove any specified sites from both surfaces.
1511
            (4) Add any specified sites to both surfaces.
1512

1513
        Returns:
1514
            (Slab): The reconstructed slab.
1515
        """
1516
        slabs = self.get_unreconstructed_slabs()
1✔
1517
        recon_slabs = []
1✔
1518

1519
        for slab in slabs:
1✔
1520
            d = get_d(slab)
1✔
1521
            top_site = sorted(slab, key=lambda site: site.frac_coords[2])[-1].coords
1✔
1522

1523
            # Remove any specified sites
1524
            if "points_to_remove" in self.reconstruction_json:
1✔
1525
                pts_to_rm = copy.deepcopy(self.reconstruction_json["points_to_remove"])
1✔
1526
                for p in pts_to_rm:
1✔
1527
                    p[2] = slab.lattice.get_fractional_coords([top_site[0], top_site[1], top_site[2] + p[2] * d])[2]
1✔
1528
                    cart_point = slab.lattice.get_cartesian_coords(p)
1✔
1529
                    dist = [site.distance_from_point(cart_point) for site in slab]
1✔
1530
                    site1 = dist.index(min(dist))
1✔
1531
                    slab.symmetrically_remove_atoms([site1])
1✔
1532

1533
            # Add any specified sites
1534
            if "points_to_add" in self.reconstruction_json:
1✔
1535
                pts_to_add = copy.deepcopy(self.reconstruction_json["points_to_add"])
1✔
1536
                for p in pts_to_add:
1✔
1537
                    p[2] = slab.lattice.get_fractional_coords([top_site[0], top_site[1], top_site[2] + p[2] * d])[2]
1✔
1538
                    slab.symmetrically_add_atom(slab[0].specie, p)
1✔
1539

1540
            slab.reconstruction = self.name
1✔
1541
            slab.recon_trans_matrix = self.trans_matrix
1✔
1542

1543
            # Get the oriented_unit_cell with the same axb area.
1544
            ouc = slab.oriented_unit_cell.copy()
1✔
1545
            ouc.make_supercell(self.trans_matrix)
1✔
1546
            slab.oriented_unit_cell = ouc
1✔
1547
            recon_slabs.append(slab)
1✔
1548

1549
        return recon_slabs
1✔
1550

1551
    def get_unreconstructed_slabs(self):
1✔
1552
        """
1553
        Generates the unreconstructed or pristine super slab.
1554
        """
1555
        slabs = []
1✔
1556
        for slab in SlabGenerator(**self.slabgen_params).get_slabs():
1✔
1557
            slab.make_supercell(self.trans_matrix)
1✔
1558
            slabs.append(slab)
1✔
1559
        return slabs
1✔
1560

1561

1562
def get_d(slab):
1✔
1563
    """
1564
    Determine the distance of space between
1565
    each layer of atoms along c
1566
    """
1567
    sorted_sites = sorted(slab, key=lambda site: site.frac_coords[2])
1✔
1568
    for i, site in enumerate(sorted_sites):
1✔
1569
        if not f"{site.frac_coords[2]:.6f}" == f"{sorted_sites[i + 1].frac_coords[2]:.6f}":
1✔
1570
            d = abs(site.frac_coords[2] - sorted_sites[i + 1].frac_coords[2])
1✔
1571
            break
1✔
1572
    return slab.lattice.get_cartesian_coords([0, 0, d])[2]
1✔
1573

1574

1575
def is_already_analyzed(miller_index: tuple, miller_list: list, symm_ops: list) -> bool:
1✔
1576
    """
1577
    Helper function to check if a given Miller index is
1578
    part of the family of indices of any index in a list
1579

1580
    Args:
1581
        miller_index (tuple): The Miller index to analyze
1582
        miller_list (list): List of Miller indices. If the given
1583
            Miller index belongs in the same family as any of the
1584
            indices in this list, return True, else return False
1585
        symm_ops (list): Symmetry operations of a
1586
            lattice, used to define family of indices
1587
    """
1588
    for op in symm_ops:
1✔
1589
        if in_coord_list(miller_list, op.operate(miller_index)):
1✔
1590
            return True
1✔
1591
    return False
1✔
1592

1593

1594
def get_symmetrically_equivalent_miller_indices(structure, miller_index, return_hkil=True):
1✔
1595
    """
1596
    Returns all symmetrically equivalent indices for a given structure. Analysis
1597
    is based on the symmetry of the reciprocal lattice of the structure.
1598

1599
    Args:
1600
        miller_index (tuple): Designates the family of Miller indices
1601
            to find. Can be hkl or hkil for hexagonal systems
1602
        return_hkil (bool): If true, return hkil form of Miller
1603
            index for hexagonal systems, otherwise return hkl
1604
    """
1605

1606
    # Change to hkl if hkil because in_coord_list only handles tuples of 3
1607
    miller_index = (miller_index[0], miller_index[1], miller_index[3]) if len(miller_index) == 4 else miller_index
1✔
1608
    mmi = max(np.abs(miller_index))
1✔
1609
    r = list(range(-mmi, mmi + 1))
1✔
1610
    r.reverse()
1✔
1611

1612
    sg = SpacegroupAnalyzer(structure)
1✔
1613
    # Get distinct hkl planes from the rhombohedral setting if trigonal
1614
    if sg.get_crystal_system() == "trigonal":
1✔
1615
        prim_structure = SpacegroupAnalyzer(structure).get_primitive_standard_structure()
×
1616
        symm_ops = prim_structure.lattice.get_recp_symmetry_operation()
×
1617
    else:
1618
        symm_ops = structure.lattice.get_recp_symmetry_operation()
1✔
1619

1620
    equivalent_millers = [miller_index]
1✔
1621
    for miller in itertools.product(r, r, r):
1✔
1622
        if miller == miller_index:
1✔
1623
            continue
1✔
1624
        if any(i != 0 for i in miller):
1✔
1625
            if is_already_analyzed(miller, equivalent_millers, symm_ops):
1✔
1626
                equivalent_millers.append(miller)
1✔
1627

1628
            # include larger Miller indices in the family of planes
1629
            if all(mmi > i for i in np.abs(miller)) and not in_coord_list(equivalent_millers, miller):
1✔
1630
                if is_already_analyzed(mmi * np.array(miller), equivalent_millers, symm_ops):
1✔
1631
                    equivalent_millers.append(miller)
1✔
1632

1633
    if return_hkil and sg.get_crystal_system() in ["trigonal", "hexagonal"]:
1✔
1634
        return [(hkl[0], hkl[1], -1 * hkl[0] - hkl[1], hkl[2]) for hkl in equivalent_millers]
1✔
1635
    return equivalent_millers
1✔
1636

1637

1638
def get_symmetrically_distinct_miller_indices(structure, max_index, return_hkil=False):
1✔
1639
    """
1640
    Returns all symmetrically distinct indices below a certain max-index for
1641
    a given structure. Analysis is based on the symmetry of the reciprocal
1642
    lattice of the structure.
1643
    Args:
1644
        structure (Structure): input structure.
1645
        max_index (int): The maximum index. For example, a max_index of 1
1646
            means that (100), (110), and (111) are returned for the cubic
1647
            structure. All other indices are equivalent to one of these.
1648
        return_hkil (bool): If true, return hkil form of Miller
1649
            index for hexagonal systems, otherwise return hkl
1650
    """
1651

1652
    r = list(range(-max_index, max_index + 1))
1✔
1653
    r.reverse()
1✔
1654

1655
    # First we get a list of all hkls for conventional (including equivalent)
1656
    conv_hkl_list = [miller for miller in itertools.product(r, r, r) if any(i != 0 for i in miller)]
1✔
1657

1658
    sg = SpacegroupAnalyzer(structure)
1✔
1659
    # Get distinct hkl planes from the rhombohedral setting if trigonal
1660
    if sg.get_crystal_system() == "trigonal":
1✔
1661
        transf = sg.get_conventional_to_primitive_transformation_matrix()
1✔
1662
        miller_list = [hkl_transformation(transf, hkl) for hkl in conv_hkl_list]
1✔
1663
        prim_structure = SpacegroupAnalyzer(structure).get_primitive_standard_structure()
1✔
1664
        symm_ops = prim_structure.lattice.get_recp_symmetry_operation()
1✔
1665
    else:
1666
        miller_list = conv_hkl_list
1✔
1667
        symm_ops = structure.lattice.get_recp_symmetry_operation()
1✔
1668

1669
    unique_millers, unique_millers_conv = [], []
1✔
1670

1671
    for i, miller in enumerate(miller_list):
1✔
1672
        d = abs(reduce(gcd, miller))
1✔
1673
        miller = tuple(int(i / d) for i in miller)
1✔
1674
        if not is_already_analyzed(miller, unique_millers, symm_ops):
1✔
1675
            if sg.get_crystal_system() == "trigonal":
1✔
1676
                # Now we find the distinct primitive hkls using
1677
                # the primitive symmetry operations and their
1678
                # corresponding hkls in the conventional setting
1679
                unique_millers.append(miller)
1✔
1680
                d = abs(reduce(gcd, conv_hkl_list[i]))
1✔
1681
                cmiller = tuple(int(i / d) for i in conv_hkl_list[i])
1✔
1682
                unique_millers_conv.append(cmiller)
1✔
1683
            else:
1684
                unique_millers.append(miller)
1✔
1685
                unique_millers_conv.append(miller)
1✔
1686

1687
    if return_hkil and sg.get_crystal_system() in ["trigonal", "hexagonal"]:
1✔
1688
        return [(hkl[0], hkl[1], -1 * hkl[0] - hkl[1], hkl[2]) for hkl in unique_millers_conv]
1✔
1689
    return unique_millers_conv
1✔
1690

1691

1692
def hkl_transformation(transf, miller_index):
1✔
1693
    """
1694
    Returns the Miller index from setting
1695
    A to B using a transformation matrix
1696
    Args:
1697
        transf (3x3 array): The transformation matrix
1698
            that transforms a lattice of A to B
1699
        miller_index ([h, k, l]): Miller index to transform to setting B
1700
    """
1701
    # Get a matrix of whole numbers (ints)
1702

1703
    def lcm(a, b):
1✔
1704
        return a * b // math.gcd(a, b)
1✔
1705

1706
    reduced_transf = reduce(lcm, [int(1 / i) for i in itertools.chain(*transf) if i != 0]) * transf
1✔
1707
    reduced_transf = reduced_transf.astype(int)
1✔
1708

1709
    # perform the transformation
1710
    t_hkl = np.dot(reduced_transf, miller_index)
1✔
1711
    d = abs(reduce(gcd, t_hkl))
1✔
1712
    t_hkl = np.array([int(i / d) for i in t_hkl])
1✔
1713

1714
    # get mostly positive oriented Miller index
1715
    if len([i for i in t_hkl if i < 0]) > 1:
1✔
1716
        t_hkl *= -1
1✔
1717

1718
    return tuple(t_hkl)
1✔
1719

1720

1721
def generate_all_slabs(
1✔
1722
    structure,
1723
    max_index,
1724
    min_slab_size,
1725
    min_vacuum_size,
1726
    bonds=None,
1727
    tol=0.1,
1728
    ftol=0.1,
1729
    max_broken_bonds=0,
1730
    lll_reduce=False,
1731
    center_slab=False,
1732
    primitive=True,
1733
    max_normal_search=None,
1734
    symmetrize=False,
1735
    repair=False,
1736
    include_reconstructions=False,
1737
    in_unit_planes=False,
1738
):
1739
    """
1740
    A function that finds all different slabs up to a certain miller index.
1741
    Slabs oriented under certain Miller indices that are equivalent to other
1742
    slabs in other Miller indices are filtered out using symmetry operations
1743
    to get rid of any repetitive slabs. For example, under symmetry operations,
1744
    CsCl has equivalent slabs in the (0,0,1), (0,1,0), and (1,0,0) direction.
1745

1746
    Args:
1747
        structure (Structure): Initial input structure. Note that to
1748
                ensure that the miller indices correspond to usual
1749
                crystallographic definitions, you should supply a conventional
1750
                unit cell structure.
1751
        max_index (int): The maximum Miller index to go up to.
1752
        min_slab_size (float): In Angstroms
1753
        min_vacuum_size (float): In Angstroms
1754
        bonds ({(specie1, specie2): max_bond_dist}: bonds are
1755
            specified as a dict of tuples: float of specie1, specie2
1756
            and the max bonding distance. For example, PO4 groups may be
1757
            defined as {("P", "O"): 3}.
1758
        tol (float): General tolerance parameter for getting primitive
1759
            cells and matching structures
1760
        ftol (float): Threshold parameter in fcluster in order to check
1761
            if two atoms are lying on the same plane. Default thresh set
1762
            to 0.1 Angstrom in the direction of the surface normal.
1763
        max_broken_bonds (int): Maximum number of allowable broken bonds
1764
            for the slab. Use this to limit # of slabs (some structures
1765
            may have a lot of slabs). Defaults to zero, which means no
1766
            defined bonds must be broken.
1767
        lll_reduce (bool): Whether to perform an LLL reduction on the
1768
            eventual structure.
1769
        center_slab (bool): Whether to center the slab in the cell with
1770
            equal vacuum spacing from the top and bottom.
1771
        primitive (bool): Whether to reduce any generated slabs to a
1772
            primitive cell (this does **not** mean the slab is generated
1773
            from a primitive cell, it simply means that after slab
1774
            generation, we attempt to find shorter lattice vectors,
1775
            which lead to less surface area and smaller cells).
1776
        max_normal_search (int): If set to a positive integer, the code will
1777
            conduct a search for a normal lattice vector that is as
1778
            perpendicular to the surface as possible by considering
1779
            multiples linear combinations of lattice vectors up to
1780
            max_normal_search. This has no bearing on surface energies,
1781
            but may be useful as a preliminary step to generating slabs
1782
            for absorption and other sizes. It is typical that this will
1783
            not be the smallest possible cell for simulation. Normality
1784
            is not guaranteed, but the oriented cell will have the c
1785
            vector as normal as possible (within the search range) to the
1786
            surface. A value of up to the max absolute Miller index is
1787
            usually sufficient.
1788
        symmetrize (bool): Whether or not to ensure the surfaces of the
1789
            slabs are equivalent.
1790
        repair (bool): Whether to repair terminations with broken bonds
1791
            or just omit them
1792
        include_reconstructions (bool): Whether to include reconstructed
1793
            slabs available in the reconstructions_archive.json file.
1794
    """
1795
    all_slabs = []
1✔
1796

1797
    for miller in get_symmetrically_distinct_miller_indices(structure, max_index):
1✔
1798
        gen = SlabGenerator(
1✔
1799
            structure,
1800
            miller,
1801
            min_slab_size,
1802
            min_vacuum_size,
1803
            lll_reduce=lll_reduce,
1804
            center_slab=center_slab,
1805
            primitive=primitive,
1806
            max_normal_search=max_normal_search,
1807
            in_unit_planes=in_unit_planes,
1808
        )
1809
        slabs = gen.get_slabs(
1✔
1810
            bonds=bonds,
1811
            tol=tol,
1812
            ftol=ftol,
1813
            symmetrize=symmetrize,
1814
            max_broken_bonds=max_broken_bonds,
1815
            repair=repair,
1816
        )
1817

1818
        if len(slabs) > 0:
1✔
1819
            logger.debug(f"{miller} has {len(slabs)} slabs... ")
1✔
1820
            all_slabs.extend(slabs)
1✔
1821

1822
    if include_reconstructions:
1✔
1823
        sg = SpacegroupAnalyzer(structure)
1✔
1824
        symbol = sg.get_space_group_symbol()
1✔
1825
        # enumerate through all posisble reconstructions in the
1826
        # archive available for this particular structure (spacegroup)
1827
        for name, instructions in reconstructions_archive.items():
1✔
1828
            if "base_reconstruction" in instructions:
1✔
1829
                instructions = reconstructions_archive[instructions["base_reconstruction"]]
1✔
1830
            if instructions["spacegroup"]["symbol"] == symbol:
1✔
1831
                # check if this reconstruction has a max index
1832
                # equal or less than the given max index
1833
                if max(instructions["miller_index"]) > max_index:
1✔
1834
                    continue
×
1835
                recon = ReconstructionGenerator(structure, min_slab_size, min_vacuum_size, name)
1✔
1836
                all_slabs.extend(recon.build_slabs())
1✔
1837

1838
    return all_slabs
1✔
1839

1840

1841
def get_slab_regions(slab, blength=3.5):
1✔
1842
    """
1843
    Function to get the ranges of the slab regions. Useful for discerning where
1844
    the slab ends and vacuum begins if the slab is not fully within the cell
1845
    Args:
1846
        slab (Structure): Structure object modelling the surface
1847
        blength (float, Ang): The bondlength between atoms. You generally
1848
            want this value to be larger than the actual bondlengths in
1849
            order to find atoms that are part of the slab
1850
    """
1851

1852
    fcoords, indices, all_indices = [], [], []
1✔
1853
    for site in slab:
1✔
1854
        # find sites with c < 0 (noncontiguous)
1855
        neighbors = slab.get_neighbors(site, blength, include_index=True, include_image=True)
1✔
1856
        for nn in neighbors:
1✔
1857
            if nn[0].frac_coords[2] < 0:
1✔
1858
                # sites are noncontiguous within cell
1859
                fcoords.append(nn[0].frac_coords[2])
1✔
1860
                indices.append(nn[-2])
1✔
1861
                if nn[-2] not in all_indices:
1✔
1862
                    all_indices.append(nn[-2])
1✔
1863

1864
    if fcoords:
1✔
1865
        # If slab is noncontiguous, locate the lowest
1866
        # site within the upper region of the slab
1867
        while fcoords:
1✔
1868
            last_fcoords = copy.copy(fcoords)
1✔
1869
            last_indices = copy.copy(indices)
1✔
1870
            site = slab[indices[fcoords.index(min(fcoords))]]
1✔
1871
            neighbors = slab.get_neighbors(site, blength, include_index=True, include_image=True)
1✔
1872
            fcoords, indices = [], []
1✔
1873
            for nn in neighbors:
1✔
1874
                if 1 > nn[0].frac_coords[2] > 0 and nn[0].frac_coords[2] < site.frac_coords[2]:
1✔
1875
                    # sites are noncontiguous within cell
1876
                    fcoords.append(nn[0].frac_coords[2])
1✔
1877
                    indices.append(nn[-2])
1✔
1878
                    if nn[-2] not in all_indices:
1✔
1879
                        all_indices.append(nn[-2])
1✔
1880

1881
        # Now locate the highest site within the lower region of the slab
1882
        upper_fcoords = []
1✔
1883
        for site in slab:
1✔
1884
            if all(nn.index not in all_indices for nn in slab.get_neighbors(site, blength)):
1✔
1885
                upper_fcoords.append(site.frac_coords[2])
1✔
1886
        coords = copy.copy(last_fcoords) if not fcoords else copy.copy(fcoords)
1✔
1887
        min_top = slab[last_indices[coords.index(min(coords))]].frac_coords[2]
1✔
1888
        ranges = [[0, max(upper_fcoords)], [min_top, 1]]
1✔
1889
    else:
1890
        # If the entire slab region is within the slab cell, just
1891
        # set the range as the highest and lowest site in the slab
1892
        sorted_sites = sorted(slab, key=lambda site: site.frac_coords[2])
1✔
1893
        ranges = [[sorted_sites[0].frac_coords[2], sorted_sites[-1].frac_coords[2]]]
1✔
1894

1895
    return ranges
1✔
1896

1897

1898
def miller_index_from_sites(lattice, coords, coords_are_cartesian=True, round_dp=4, verbose=True):
1✔
1899
    """
1900
    Get the Miller index of a plane from a list of site coordinates.
1901

1902
    A minimum of 3 sets of coordinates are required. If more than 3 sets of
1903
    coordinates are given, the best plane that minimises the distance to all
1904
    points will be calculated.
1905

1906
    Args:
1907
        lattice (list or Lattice): A 3x3 lattice matrix or `Lattice` object (for
1908
            example obtained from Structure.lattice).
1909
        coords (iterable): A list or numpy array of coordinates. Can be
1910
            Cartesian or fractional coordinates. If more than three sets of
1911
            coordinates are provided, the best plane that minimises the
1912
            distance to all sites will be calculated.
1913
        coords_are_cartesian (bool, optional): Whether the coordinates are
1914
            in Cartesian space. If using fractional coordinates set to False.
1915
        round_dp (int, optional): The number of decimal places to round the
1916
            miller index to.
1917
        verbose (bool, optional): Whether to print warnings.
1918

1919
    Returns:
1920
        (tuple): The Miller index.
1921
    """
1922
    if not isinstance(lattice, Lattice):
1✔
1923
        lattice = Lattice(lattice)
1✔
1924

1925
    return lattice.get_miller_index_from_coords(
1✔
1926
        coords,
1927
        coords_are_cartesian=coords_are_cartesian,
1928
        round_dp=round_dp,
1929
        verbose=verbose,
1930
    )
1931

1932

1933
def center_slab(slab):
1✔
1934
    """
1935
    The goal here is to ensure the center of the slab region
1936
        is centered close to c=0.5. This makes it easier to
1937
        find the surface sites and apply operations like doping.
1938

1939
    There are three cases where the slab in not centered:
1940

1941
    1. The slab region is completely between two vacuums in the
1942
    box but not necessarily centered. We simply shift the
1943
    slab by the difference in its center of mass and 0.5
1944
    along the c direction.
1945

1946
    2. The slab completely spills outside the box from the bottom
1947
    and into the top. This makes it incredibly difficult to
1948
    locate surface sites. We iterate through all sites that
1949
    spill over (z>c) and shift all sites such that this specific
1950
    site is now on the other side. Repeat for all sites with z>c.
1951

1952
    3. This is a simpler case of scenario 2. Either the top or bottom
1953
    slab sites are at c=0 or c=1. Treat as scenario 2.
1954

1955
    Args:
1956
        slab (Slab): Slab structure to center
1957
    Returns:
1958
        Returns a centered slab structure
1959
    """
1960

1961
    # get a reasonable r cutoff to sample neighbors
1962
    bdists = sorted(nn[1] for nn in slab.get_neighbors(slab[0], 10) if nn[1] > 0)
×
1963
    r = bdists[0] * 3
×
1964

1965
    all_indices = [i for i, site in enumerate(slab)]
×
1966

1967
    # check if structure is case 2 or 3, shift all the
1968
    # sites up to the other side until it is case 1
1969
    for site in slab:
×
1970
        if any(nn[1] > slab.lattice.c for nn in slab.get_neighbors(site, r)):
×
1971
            shift = 1 - site.frac_coords[2] + 0.05
×
1972
            slab.translate_sites(all_indices, [0, 0, shift])
×
1973

1974
    # now the slab is case 1, shift the center of mass of the slab to 0.5
1975
    weights = [s.species.weight for s in slab]
×
1976
    center_of_mass = np.average(slab.frac_coords, weights=weights, axis=0)
×
1977
    shift = 0.5 - center_of_mass[2]
×
1978
    slab.translate_sites(all_indices, [0, 0, shift])
×
1979

1980
    return slab
×
1981

1982

1983
def _reduce_vector(vector):
1✔
1984
    # small function to reduce vectors
1985

1986
    d = abs(reduce(gcd, vector))
1✔
1987
    vector = tuple(int(i / d) for i in vector)
1✔
1988

1989
    return vector
1✔
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2025 Coveralls, Inc