• 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

80.76
/pymatgen/analysis/adsorption.py
1
# Copyright (c) Pymatgen Development Team.
2
# Distributed under the terms of the MIT License.
3

4
"""This module provides classes used to enumerate surface sites and to find
1✔
5
adsorption sites on slabs."""
6

7
from __future__ import annotations
1✔
8

9
import itertools
1✔
10
import os
1✔
11

12
import numpy as np
1✔
13
from matplotlib import patches
1✔
14
from matplotlib.path import Path
1✔
15
from monty.serialization import loadfn
1✔
16
from scipy.spatial import Delaunay
1✔
17

18
from pymatgen import vis
1✔
19
from pymatgen.analysis.local_env import VoronoiNN
1✔
20
from pymatgen.analysis.structure_matcher import StructureMatcher
1✔
21
from pymatgen.core.operations import SymmOp
1✔
22
from pymatgen.core.structure import Structure
1✔
23
from pymatgen.core.surface import generate_all_slabs
1✔
24
from pymatgen.symmetry.analyzer import SpacegroupAnalyzer
1✔
25
from pymatgen.util.coord import in_coord_list_pbc
1✔
26

27
__author__ = "Joseph Montoya"
1✔
28
__copyright__ = "Copyright 2016, The Materials Project"
1✔
29
__version__ = "0.1"
1✔
30
__maintainer__ = "Joseph Montoya"
1✔
31
__credits__ = "Richard Tran"
1✔
32
__email__ = "montoyjh@lbl.gov"
1✔
33
__status__ = "Development"
1✔
34
__date__ = "December 2, 2015"
1✔
35

36

37
class AdsorbateSiteFinder:
1✔
38
    """This class finds adsorbate sites on slabs and generates adsorbate
39
    structures according to user-defined criteria.
40

41
    The algorithm for finding sites is essentially as follows:
42
        1. Determine "surface sites" by finding those within
43
            a height threshold along the miller index of the
44
            highest site
45
        2. Create a network of surface sites using the Delaunay
46
            triangulation of the surface sites
47
        3. Assign on-top, bridge, and hollow adsorption sites
48
            at the nodes, edges, and face centers of the Del.
49
            Triangulation
50
        4. Generate structures from a molecule positioned at
51
            these sites
52
    """
53

54
    def __init__(self, slab, selective_dynamics=False, height=0.9, mi_vec=None):
1✔
55
        """Create an AdsorbateSiteFinder object.
56

57
        Args:
58
            slab (Slab): slab object for which to find adsorbate sites
59
            selective_dynamics (bool): flag for whether to assign
60
                non-surface sites as fixed for selective dynamics
61
            height (float): height criteria for selection of surface sites
62
            mi_vec (3-D array-like): vector corresponding to the vector
63
                concurrent with the miller index, this enables use with
64
                slabs that have been reoriented, but the miller vector
65
                must be supplied manually
66
        """
67
        # get surface normal from miller index
68
        if mi_vec:
1✔
69
            self.mvec = mi_vec
×
70
        else:
71
            self.mvec = get_mi_vec(slab)
1✔
72
        slab = self.assign_site_properties(slab, height)
1✔
73
        if selective_dynamics:
1✔
74
            slab = self.assign_selective_dynamics(slab)
×
75
        self.slab = slab
1✔
76

77
    @classmethod
1✔
78
    def from_bulk_and_miller(
1✔
79
        cls,
80
        structure,
81
        miller_index,
82
        min_slab_size=8.0,
83
        min_vacuum_size=10.0,
84
        max_normal_search=None,
85
        center_slab=True,
86
        selective_dynamics=False,
87
        undercoord_threshold=0.09,
88
    ):
89
        """This method constructs the adsorbate site finder from a bulk
90
        structure and a miller index, which allows the surface sites to be
91
        determined from the difference in bulk and slab coordination, as
92
        opposed to the height threshold.
93

94
        Args:
95
            structure (Structure): structure from which slab
96
                input to the ASF is constructed
97
            miller_index (3-tuple or list): miller index to be used
98
            min_slab_size (float): min slab size for slab generation
99
            min_vacuum_size (float): min vacuum size for slab generation
100
            max_normal_search (int): max normal search for slab generation
101
            center_slab (bool): whether to center slab in slab generation
102
            selective dynamics (bool): whether to assign surface sites
103
                to selective dynamics
104
            undercoord_threshold (float): threshold of "undercoordation"
105
                to use for the assignment of surface sites. Default is
106
                0.1, for which surface sites will be designated if they
107
                are 10% less coordinated than their bulk counterpart
108
        """
109
        # TODO: for some reason this works poorly with primitive cells
110
        #       may want to switch the coordination algorithm eventually
111
        vnn_bulk = VoronoiNN(tol=0.05)
1✔
112
        bulk_coords = [len(vnn_bulk.get_nn(structure, n)) for n in range(len(structure))]
1✔
113
        struct = structure.copy(site_properties={"bulk_coordinations": bulk_coords})
1✔
114
        slabs = generate_all_slabs(
1✔
115
            struct,
116
            max_index=max(miller_index),
117
            min_slab_size=min_slab_size,
118
            min_vacuum_size=min_vacuum_size,
119
            max_normal_search=max_normal_search,
120
            center_slab=center_slab,
121
        )
122

123
        slab_dict = {slab.miller_index: slab for slab in slabs}
1✔
124

125
        if miller_index not in slab_dict:
1✔
126
            raise ValueError("Miller index not in slab dict")
×
127

128
        this_slab = slab_dict[miller_index]
1✔
129

130
        vnn_surface = VoronoiNN(tol=0.05, allow_pathological=True)
1✔
131

132
        surf_props, undercoords = [], []
1✔
133
        this_mi_vec = get_mi_vec(this_slab)
1✔
134
        mi_mags = [np.dot(this_mi_vec, site.coords) for site in this_slab]
1✔
135
        average_mi_mag = np.average(mi_mags)
1✔
136
        for n, site in enumerate(this_slab):
1✔
137
            bulk_coord = this_slab.site_properties["bulk_coordinations"][n]
1✔
138
            slab_coord = len(vnn_surface.get_nn(this_slab, n))
1✔
139
            mi_mag = np.dot(this_mi_vec, site.coords)
1✔
140
            undercoord = (bulk_coord - slab_coord) / bulk_coord
1✔
141
            undercoords += [undercoord]
1✔
142
            if undercoord > undercoord_threshold and mi_mag > average_mi_mag:
1✔
143
                surf_props += ["surface"]
1✔
144
            else:
145
                surf_props += ["subsurface"]
1✔
146
        new_site_properties = {
1✔
147
            "surface_properties": surf_props,
148
            "undercoords": undercoords,
149
        }
150
        new_slab = this_slab.copy(site_properties=new_site_properties)
1✔
151
        return cls(new_slab, selective_dynamics)
1✔
152

153
    def find_surface_sites_by_height(self, slab, height=0.9, xy_tol=0.05):
1✔
154
        """This method finds surface sites by determining which sites are
155
        within a threshold value in height from the topmost site in a list of
156
        sites.
157

158
        Args:
159
            site_list (list): list of sites from which to select surface sites
160
            height (float): threshold in angstroms of distance from topmost
161
                site in slab along the slab c-vector to include in surface
162
                site determination
163
            xy_tol (float): if supplied, will remove any sites which are
164
                within a certain distance in the miller plane.
165

166
        Returns:
167
            list of sites selected to be within a threshold of the highest
168
        """
169
        # Get projection of coordinates along the miller index
170
        m_projs = np.array([np.dot(site.coords, self.mvec) for site in slab.sites])
1✔
171

172
        # Mask based on window threshold along the miller index.
173
        mask = (m_projs - np.amax(m_projs)) >= -height
1✔
174
        surf_sites = [slab.sites[n] for n in np.where(mask)[0]]
1✔
175
        if xy_tol:
1✔
176
            # sort surface sites by height
177
            surf_sites = [s for (h, s) in zip(m_projs[mask], surf_sites)]
1✔
178
            surf_sites.reverse()
1✔
179
            unique_sites, unique_perp_fracs = [], []
1✔
180
            for site in surf_sites:
1✔
181
                this_perp = site.coords - np.dot(site.coords, self.mvec)
1✔
182
                this_perp_frac = slab.lattice.get_fractional_coords(this_perp)
1✔
183
                if not in_coord_list_pbc(unique_perp_fracs, this_perp_frac):
1✔
184
                    unique_sites.append(site)
1✔
185
                    unique_perp_fracs.append(this_perp_frac)
1✔
186
            surf_sites = unique_sites
1✔
187

188
        return surf_sites
1✔
189

190
    def assign_site_properties(self, slab, height=0.9):
1✔
191
        """Assigns site properties."""
192
        if "surface_properties" in slab.site_properties:
1✔
193
            return slab
1✔
194

195
        surf_sites = self.find_surface_sites_by_height(slab, height)
1✔
196
        surf_props = ["surface" if site in surf_sites else "subsurface" for site in slab.sites]
1✔
197
        return slab.copy(site_properties={"surface_properties": surf_props})
1✔
198

199
    def get_extended_surface_mesh(self, repeat=(5, 5, 1)):
1✔
200
        """Gets an extended surface mesh for to use for adsorption site finding
201
        by constructing supercell of surface sites.
202

203
        Args:
204
            repeat (3-tuple): repeat for getting extended surface mesh
205
        """
206
        surf_str = Structure.from_sites(self.surface_sites)
1✔
207
        surf_str.make_supercell(repeat)
1✔
208
        return surf_str
1✔
209

210
    @property
1✔
211
    def surface_sites(self):
1✔
212
        """convenience method to return a list of surface sites."""
213
        return [site for site in self.slab.sites if site.properties["surface_properties"] == "surface"]
1✔
214

215
    def subsurface_sites(self):
1✔
216
        """convenience method to return list of subsurface sites."""
217
        return [site for site in self.slab.sites if site.properties["surface_properties"] == "subsurface"]
1✔
218

219
    def find_adsorption_sites(
1✔
220
        self,
221
        distance=2.0,
222
        put_inside=True,
223
        symm_reduce=1e-2,
224
        near_reduce=1e-2,
225
        positions=("ontop", "bridge", "hollow"),
226
        no_obtuse_hollow=True,
227
    ):
228
        """Finds surface sites according to the above algorithm. Returns a list
229
        of corresponding Cartesian coordinates.
230

231
        Args:
232
            distance (float): distance from the coordinating ensemble
233
                of atoms along the miller index for the site (i. e.
234
                the distance from the slab itself)
235
            put_inside (bool): whether to put the site inside the cell
236
            symm_reduce (float): symm reduction threshold
237
            near_reduce (float): near reduction threshold
238
            positions (list): which positions to include in the site finding
239
                "ontop": sites on top of surface sites
240
                "bridge": sites at edges between surface sites in Delaunay
241
                    triangulation of surface sites in the miller plane
242
                "hollow": sites at centers of Delaunay triangulation faces
243
                "subsurface": subsurface positions projected into miller plane
244
            no_obtuse_hollow (bool): flag to indicate whether to include
245
                obtuse triangular ensembles in hollow sites
246
        """
247
        ads_sites = {k: [] for k in positions}
1✔
248
        if "ontop" in positions:
1✔
249
            ads_sites["ontop"] = [s.coords for s in self.surface_sites]
1✔
250
        if "subsurface" in positions:
1✔
251
            # Get highest site
252
            ref = self.slab.sites[np.argmax(self.slab.cart_coords[:, 2])]
1✔
253
            # Project diff between highest site and subs site into miller
254
            ss_sites = [
1✔
255
                self.mvec * np.dot(ref.coords - s.coords, self.mvec) + s.coords for s in self.subsurface_sites()
256
            ]
257
            ads_sites["subsurface"] = ss_sites
1✔
258
        if "bridge" in positions or "hollow" in positions:
1✔
259
            mesh = self.get_extended_surface_mesh()
1✔
260
            sop = get_rot(self.slab)
1✔
261
            dt = Delaunay([sop.operate(m.coords)[:2] for m in mesh])
1✔
262
            # TODO: refactor below to properly account for >3-fold
263
            for v in dt.simplices:
1✔
264
                if -1 not in v:
1✔
265
                    dots = []
1✔
266
                    for i_corner, i_opp in zip(range(3), ((1, 2), (0, 2), (0, 1))):
1✔
267
                        corner, opp = v[i_corner], [v[o] for o in i_opp]
1✔
268
                        vecs = [mesh[d].coords - mesh[corner].coords for d in opp]
1✔
269
                        vecs = [vec / np.linalg.norm(vec) for vec in vecs]
1✔
270
                        dots.append(np.dot(*vecs))
1✔
271
                        # Add bridge sites at midpoints of edges of D. Tri
272
                        if "bridge" in positions:
1✔
273
                            ads_sites["bridge"].append(self.ensemble_center(mesh, opp))
1✔
274
                    # Prevent addition of hollow sites in obtuse triangles
275
                    obtuse = no_obtuse_hollow and (np.array(dots) < 1e-5).any()
1✔
276
                    # Add hollow sites at centers of D. Tri faces
277
                    if "hollow" in positions and not obtuse:
1✔
278
                        ads_sites["hollow"].append(self.ensemble_center(mesh, v))
1✔
279
        for key, sites in ads_sites.items():
1✔
280
            # Pare off outer sites for bridge/hollow
281
            if key in ["bridge", "hollow"]:
1✔
282
                frac_coords = [self.slab.lattice.get_fractional_coords(ads_site) for ads_site in sites]
1✔
283
                frac_coords = [
1✔
284
                    frac_coord
285
                    for frac_coord in frac_coords
286
                    if (frac_coord[0] > 1 and frac_coord[0] < 4 and frac_coord[1] > 1 and frac_coord[1] < 4)
287
                ]
288
                sites = [self.slab.lattice.get_cartesian_coords(frac_coord) for frac_coord in frac_coords]
1✔
289
            if near_reduce:
1✔
290
                sites = self.near_reduce(sites, threshold=near_reduce)
1✔
291
            if put_inside:
1✔
292
                sites = [put_coord_inside(self.slab.lattice, coord) for coord in sites]
1✔
293
            if symm_reduce:
1✔
294
                sites = self.symm_reduce(sites, threshold=symm_reduce)
1✔
295
            sites = [site + distance * self.mvec for site in sites]
1✔
296

297
            ads_sites[key] = sites
1✔
298
        ads_sites["all"] = sum(ads_sites.values(), [])
1✔
299
        return ads_sites
1✔
300

301
    def symm_reduce(self, coords_set, threshold=1e-6):
1✔
302
        """Reduces the set of adsorbate sites by finding removing symmetrically
303
        equivalent duplicates.
304

305
        Args:
306
            coords_set: coordinate set in Cartesian coordinates
307
            threshold: tolerance for distance equivalence, used
308
                as input to in_coord_list_pbc for dupl. checking
309
        """
310
        surf_sg = SpacegroupAnalyzer(self.slab, 0.1)
1✔
311
        symm_ops = surf_sg.get_symmetry_operations()
1✔
312
        unique_coords = []
1✔
313
        # Convert to fractional
314
        coords_set = [self.slab.lattice.get_fractional_coords(coords) for coords in coords_set]
1✔
315
        for coords in coords_set:
1✔
316
            incoord = False
1✔
317
            for op in symm_ops:
1✔
318
                if in_coord_list_pbc(unique_coords, op.operate(coords), atol=threshold):
1✔
319
                    incoord = True
1✔
320
                    break
1✔
321
            if not incoord:
1✔
322
                unique_coords += [coords]
1✔
323
        # convert back to cartesian
324
        return [self.slab.lattice.get_cartesian_coords(coords) for coords in unique_coords]
1✔
325

326
    def near_reduce(self, coords_set, threshold=1e-4):
1✔
327
        """Prunes coordinate set for coordinates that are within threshold.
328

329
        Args:
330
            coords_set (Nx3 array-like): list or array of coordinates
331
            threshold (float): threshold value for distance
332
        """
333
        unique_coords = []
1✔
334
        coords_set = [self.slab.lattice.get_fractional_coords(coords) for coords in coords_set]
1✔
335
        for coord in coords_set:
1✔
336
            if not in_coord_list_pbc(unique_coords, coord, threshold):
1✔
337
                unique_coords += [coord]
1✔
338
        return [self.slab.lattice.get_cartesian_coords(coords) for coords in unique_coords]
1✔
339

340
    @classmethod
1✔
341
    def ensemble_center(cls, site_list, indices, cartesian=True):
1✔
342
        """Finds the center of an ensemble of sites selected from a list of
343
        sites. Helper method for the find_adsorption_sites algorithm.
344

345
        Args:
346
            site_list (list of sites): list of sites
347
            indices (list of ints): list of ints from which to select
348
                sites from site list
349
            cartesian (bool): whether to get average fractional or
350
                Cartesian coordinate
351
        """
352
        if cartesian:
1✔
353
            return np.average([site_list[i].coords for i in indices], axis=0)
1✔
354

355
        return np.average([site_list[i].frac_coords for i in indices], axis=0)
×
356

357
    def add_adsorbate(self, molecule, ads_coord, repeat=None, translate=True, reorient=True):
1✔
358
        """Adds an adsorbate at a particular coordinate. Adsorbate represented
359
        by a Molecule object and is translated to (0, 0, 0) if translate is
360
        True, or positioned relative to the input adsorbate coordinate if
361
        translate is False.
362

363
        Args:
364
            molecule (Molecule): molecule object representing the adsorbate
365
            ads_coord (array): coordinate of adsorbate position
366
            repeat (3-tuple or list): input for making a supercell of slab
367
                prior to placing the adsorbate
368
            translate (bool): flag on whether to translate the molecule so
369
                that its CoM is at the origin prior to adding it to the surface
370
            reorient (bool): flag on whether to reorient the molecule to
371
                have its z-axis concurrent with miller index
372
        """
373
        molecule = molecule.copy()
1✔
374
        if translate:
1✔
375
            # Translate the molecule so that the center of mass of the atoms
376
            # that have the most negative z coordinate is at (0, 0, 0)
377
            front_atoms = molecule.copy()
1✔
378
            front_atoms._sites = [s for s in molecule.sites if s.coords[2] == min(s.coords[2] for s in molecule.sites)]
1✔
379
            x, y, z = front_atoms.center_of_mass
1✔
380
            molecule.translate_sites(vector=[-x, -y, -z])
1✔
381
        if reorient:
1✔
382
            # Reorient the molecule along slab m_index
383
            sop = get_rot(self.slab)
1✔
384
            molecule.apply_operation(sop.inverse)
1✔
385
        struct = self.slab.copy()
1✔
386
        if repeat:
1✔
387
            struct.make_supercell(repeat)
1✔
388
        if "surface_properties" in struct.site_properties:
1✔
389
            molecule.add_site_property("surface_properties", ["adsorbate"] * molecule.num_sites)
1✔
390
        if "selective_dynamics" in struct.site_properties:
1✔
391
            molecule.add_site_property("selective_dynamics", [[True, True, True]] * molecule.num_sites)
×
392
        for site in molecule:
1✔
393
            struct.append(
1✔
394
                site.specie,
395
                ads_coord + site.coords,
396
                coords_are_cartesian=True,
397
                properties=site.properties,
398
            )
399
        return struct
1✔
400

401
    @classmethod
1✔
402
    def assign_selective_dynamics(cls, slab):
1✔
403
        """Helper function to assign selective dynamics site_properties based
404
        on surface, subsurface site properties.
405

406
        Args:
407
            slab (Slab): slab for which to assign selective dynamics
408
        """
409
        sd_list = []
×
410
        sd_list = [
×
411
            [False, False, False] if site.properties["surface_properties"] == "subsurface" else [True, True, True]
412
            for site in slab.sites
413
        ]
414
        new_sp = slab.site_properties
×
415
        new_sp["selective_dynamics"] = sd_list
×
416
        return slab.copy(site_properties=new_sp)
×
417

418
    def generate_adsorption_structures(
1✔
419
        self,
420
        molecule,
421
        repeat=None,
422
        min_lw=5.0,
423
        translate=True,
424
        reorient=True,
425
        find_args=None,
426
    ):
427
        """Function that generates all adsorption structures for a given
428
        molecular adsorbate. Can take repeat argument or minimum length/width
429
        of precursor slab as an input.
430

431
        Args:
432
            molecule (Molecule): molecule corresponding to adsorbate
433
            repeat (3-tuple or list): repeat argument for supercell generation
434
            min_lw (float): minimum length and width of the slab, only used
435
                if repeat is None
436
            translate (bool): flag on whether to translate the molecule so
437
                that its CoM is at the origin prior to adding it to the surface
438
            reorient (bool): flag on whether or not to reorient adsorbate
439
                along the miller index
440
            find_args (dict): dictionary of arguments to be passed to the
441
                call to self.find_adsorption_sites, e.g. {"distance":2.0}
442
        """
443
        if repeat is None:
1✔
444
            xrep = np.ceil(min_lw / np.linalg.norm(self.slab.lattice.matrix[0]))
1✔
445
            yrep = np.ceil(min_lw / np.linalg.norm(self.slab.lattice.matrix[1]))
1✔
446
            repeat = [xrep, yrep, 1]
1✔
447
        structs = []
1✔
448

449
        find_args = find_args or {}
1✔
450
        for coords in self.find_adsorption_sites(**find_args)["all"]:
1✔
451
            structs.append(
1✔
452
                self.add_adsorbate(
453
                    molecule,
454
                    coords,
455
                    repeat=repeat,
456
                    translate=translate,
457
                    reorient=reorient,
458
                )
459
            )
460
        return structs
1✔
461

462
    def adsorb_both_surfaces(
1✔
463
        self,
464
        molecule,
465
        repeat=None,
466
        min_lw=5.0,
467
        translate=True,
468
        reorient=True,
469
        find_args=None,
470
    ):
471
        """Function that generates all adsorption structures for a given
472
        molecular adsorbate on both surfaces of a slab. This is useful for
473
        calculating surface energy where both surfaces need to be equivalent or
474
        if we want to calculate nonpolar systems.
475

476
        Args:
477
            molecule (Molecule): molecule corresponding to adsorbate
478
            repeat (3-tuple or list): repeat argument for supercell generation
479
            min_lw (float): minimum length and width of the slab, only used
480
                if repeat is None
481
            reorient (bool): flag on whether or not to reorient adsorbate
482
                along the miller index
483
            find_args (dict): dictionary of arguments to be passed to the
484
                call to self.find_adsorption_sites, e.g. {"distance":2.0}
485
        """
486
        # Get the adsorbed surfaces first
487
        find_args = find_args or {}
1✔
488
        ad_slabss = self.generate_adsorption_structures(
1✔
489
            molecule,
490
            repeat=repeat,
491
            min_lw=min_lw,
492
            translate=translate,
493
            reorient=reorient,
494
            find_args=find_args,
495
        )
496

497
        new_ad_slabss = []
1✔
498
        for ad_slabs in ad_slabss:
1✔
499
            # Find the adsorbate sites and indices in each slab
500
            _, adsorbates, indices = False, [], []
1✔
501
            for i, site in enumerate(ad_slabs.sites):
1✔
502
                if site.surface_properties == "adsorbate":
1✔
503
                    adsorbates.append(site)
1✔
504
                    indices.append(i)
1✔
505

506
            # Start with the clean slab
507
            ad_slabs.remove_sites(indices)
1✔
508
            slab = ad_slabs.copy()
1✔
509

510
            # For each site, we add it back to the slab along with a
511
            # symmetrically equivalent position on the other side of
512
            # the slab using symmetry operations
513
            for adsorbate in adsorbates:
1✔
514
                p2 = ad_slabs.get_symmetric_site(adsorbate.frac_coords)
1✔
515
                slab.append(adsorbate.specie, p2, properties={"surface_properties": "adsorbate"})
1✔
516
                slab.append(
1✔
517
                    adsorbate.specie,
518
                    adsorbate.frac_coords,
519
                    properties={"surface_properties": "adsorbate"},
520
                )
521
            new_ad_slabss.append(slab)
1✔
522

523
        return new_ad_slabss
1✔
524

525
    def generate_substitution_structures(
1✔
526
        self,
527
        atom,
528
        target_species=None,
529
        sub_both_sides=False,
530
        range_tol=1e-2,
531
        dist_from_surf=0,
532
    ):
533
        """Function that performs substitution-type doping on the surface and
534
        returns all possible configurations where one dopant is substituted per
535
        surface. Can substitute one surface or both.
536

537
        Args:
538
            atom (str): atom corresponding to substitutional dopant
539
            sub_both_sides (bool): If true, substitute an equivalent
540
                site on the other surface
541
            target_species (list): List of specific species to substitute
542
            range_tol (float): Find viable substitution sites at a specific
543
                distance from the surface +- this tolerance
544
            dist_from_surf (float): Distance from the surface to find viable
545
                substitution sites, defaults to 0 to substitute at the surface
546
        """
547
        target_species = target_species or []
1✔
548

549
        # Get symmetrized structure in case we want to substitute both sides
550
        sym_slab = SpacegroupAnalyzer(self.slab).get_symmetrized_structure()
1✔
551

552
        # Define a function for substituting a site
553
        def substitute(site, i):
1✔
554
            slab = self.slab.copy()
1✔
555
            props = self.slab.site_properties
1✔
556
            if sub_both_sides:
1✔
557
                # Find an equivalent site on the other surface
558
                eq_indices = [indices for indices in sym_slab.equivalent_indices if i in indices][0]
1✔
559
                for ii in eq_indices:
1✔
560
                    if f"{sym_slab[ii].frac_coords[2]:.6f}" != f"{site.frac_coords[2]:.6f}":
1✔
561
                        props["surface_properties"][ii] = "substitute"
1✔
562
                        slab.replace(ii, atom)
1✔
563
                        break
1✔
564

565
            props["surface_properties"][i] = "substitute"
1✔
566
            slab.replace(i, atom)
1✔
567
            slab.add_site_property("surface_properties", props["surface_properties"])
1✔
568
            return slab
1✔
569

570
        # Get all possible substitution sites
571
        substituted_slabs = []
1✔
572
        # Sort sites so that we can define a range relative to the position of the
573
        # surface atoms, i.e. search for sites above (below) the bottom (top) surface
574
        sorted_sites = sorted(sym_slab, key=lambda site: site.frac_coords[2])
1✔
575
        if sorted_sites[0].surface_properties == "surface":
1✔
576
            d = sorted_sites[0].frac_coords[2] + dist_from_surf
×
577
        else:
578
            d = sorted_sites[-1].frac_coords[2] - dist_from_surf
1✔
579

580
        for i, site in enumerate(sym_slab):
1✔
581
            if d - range_tol < site.frac_coords[2] < d + range_tol:
1✔
582
                if target_species and site.species_string in target_species:
1✔
583
                    substituted_slabs.append(substitute(site, i))
1✔
584
                elif not target_species:
1✔
585
                    substituted_slabs.append(substitute(site, i))
1✔
586

587
        matcher = StructureMatcher()
1✔
588
        return [s[0] for s in matcher.group_structures(substituted_slabs)]
1✔
589

590

591
def get_mi_vec(slab):
1✔
592
    """Convenience function which returns the unit vector aligned with the
593
    miller index."""
594
    mvec = np.cross(slab.lattice.matrix[0], slab.lattice.matrix[1])
1✔
595
    return mvec / np.linalg.norm(mvec)
1✔
596

597

598
def get_rot(slab):
1✔
599
    """Gets the transformation to rotate the z axis into the miller index."""
600
    new_z = get_mi_vec(slab)
1✔
601
    a, b, c = slab.lattice.matrix
1✔
602
    new_x = a / np.linalg.norm(a)
1✔
603
    new_y = np.cross(new_z, new_x)
1✔
604
    x, y, z = np.eye(3)
1✔
605
    rot_matrix = np.array([np.dot(*el) for el in itertools.product([x, y, z], [new_x, new_y, new_z])]).reshape(3, 3)
1✔
606
    rot_matrix = np.transpose(rot_matrix)
1✔
607
    sop = SymmOp.from_rotation_and_translation(rot_matrix)
1✔
608
    return sop
1✔
609

610

611
def put_coord_inside(lattice, cart_coordinate):
1✔
612
    """converts a Cartesian coordinate such that it is inside the unit cell."""
613
    fc = lattice.get_fractional_coords(cart_coordinate)
1✔
614
    return lattice.get_cartesian_coords([c - np.floor(c) for c in fc])
1✔
615

616

617
def reorient_z(structure):
1✔
618
    """reorients a structure such that the z axis is concurrent with the normal
619
    to the A-B plane."""
620
    struct = structure.copy()
1✔
621
    sop = get_rot(struct)
1✔
622
    struct.apply_operation(sop)
1✔
623
    return struct
1✔
624

625

626
# Get color dictionary
627
colors = loadfn(os.path.join(os.path.dirname(vis.__file__), "ElementColorSchemes.yaml"))
1✔
628
color_dict = {el: [j / 256.001 for j in colors["Jmol"][el]] for el in colors["Jmol"]}
1✔
629

630

631
def plot_slab(
1✔
632
    slab,
633
    ax,
634
    scale=0.8,
635
    repeat=5,
636
    window=1.5,
637
    draw_unit_cell=True,
638
    decay=0.2,
639
    adsorption_sites=True,
640
    inverse=False,
641
):
642
    """Function that helps visualize the slab in a 2-D plot, for convenient
643
    viewing of output of AdsorbateSiteFinder.
644

645
    Args:
646
        slab (slab): Slab object to be visualized
647
        ax (axes): matplotlib axes with which to visualize
648
        scale (float): radius scaling for sites
649
        repeat (int): number of repeating unit cells to visualize
650
        window (float): window for setting the axes limits, is essentially
651
            a fraction of the unit cell limits
652
        draw_unit_cell (bool): flag indicating whether or not to draw cell
653
        decay (float): how the alpha-value decays along the z-axis
654
        inverse (bool): invert z axis to plot opposite surface
655
    """
656
    orig_slab = slab.copy()
×
657
    slab = reorient_z(slab)
×
658
    orig_cell = slab.lattice.matrix.copy()
×
659
    if repeat:
×
660
        slab.make_supercell([repeat, repeat, 1])
×
661
    coords = np.array(sorted(slab.cart_coords, key=lambda x: x[2]))
×
662
    sites = sorted(slab.sites, key=lambda x: x.coords[2])
×
663
    alphas = 1 - decay * (np.max(coords[:, 2]) - coords[:, 2])
×
664
    alphas = alphas.clip(min=0)
×
665
    corner = [0, 0, slab.lattice.get_fractional_coords(coords[-1])[-1]]
×
666
    corner = slab.lattice.get_cartesian_coords(corner)[:2]
×
667
    verts = orig_cell[:2, :2]
×
668
    lattsum = verts[0] + verts[1]
×
669
    # inverse coords, sites, alphas, to show other side of slab
670
    if inverse:
×
671
        alphas = np.array(reversed(alphas))
×
672
        sites = list(reversed(sites))
×
673
        coords = np.array(reversed(coords))
×
674
    # Draw circles at sites and stack them accordingly
675
    for n, coord in enumerate(coords):
×
676
        r = sites[n].species.elements[0].atomic_radius * scale
×
677
        ax.add_patch(patches.Circle(coord[:2] - lattsum * (repeat // 2), r, color="w", zorder=2 * n))
×
678
        color = color_dict[sites[n].species.elements[0].symbol]
×
679
        ax.add_patch(
×
680
            patches.Circle(
681
                coord[:2] - lattsum * (repeat // 2),
682
                r,
683
                facecolor=color,
684
                alpha=alphas[n],
685
                edgecolor="k",
686
                lw=0.3,
687
                zorder=2 * n + 1,
688
            )
689
        )
690
    # Adsorption sites
691
    if adsorption_sites:
×
692
        asf = AdsorbateSiteFinder(orig_slab)
×
693
        if inverse:
×
694
            inverse_slab = orig_slab.copy()
×
695
            inverse_slab.make_supercell([1, 1, -1])
×
696
            asf = AdsorbateSiteFinder(inverse_slab)
×
697
        ads_sites = asf.find_adsorption_sites()["all"]
×
698
        sop = get_rot(orig_slab)
×
699
        ads_sites = [sop.operate(ads_site)[:2].tolist() for ads_site in ads_sites]
×
700
        ax.plot(*zip(*ads_sites), color="k", marker="x", markersize=10, mew=1, linestyle="", zorder=10000)
×
701
    # Draw unit cell
702
    if draw_unit_cell:
×
703
        verts = np.insert(verts, 1, lattsum, axis=0).tolist()
×
704
        verts += [[0.0, 0.0]]
×
705
        verts = [[0.0, 0.0]] + verts
×
706
        codes = [Path.MOVETO, Path.LINETO, Path.LINETO, Path.LINETO, Path.CLOSEPOLY]
×
707
        verts = [(np.array(vert) + corner).tolist() for vert in verts]
×
708
        path = Path(verts, codes)
×
709
        patch = patches.PathPatch(path, facecolor="none", lw=2, alpha=0.5, zorder=2 * n + 2)
×
710
        ax.add_patch(patch)
×
711
    ax.set_aspect("equal")
×
712
    center = corner + lattsum / 2.0
×
713
    extent = np.max(lattsum)
×
714
    lim_array = [center - extent * window, center + extent * window]
×
715
    x_lim = [ele[0] for ele in lim_array]
×
716
    y_lim = [ele[1] for ele in lim_array]
×
717
    ax.set_xlim(x_lim)
×
718
    ax.set_ylim(y_lim)
×
719
    return ax
×
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